Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
AnasaziTraceMinBase.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 
10 // TODO: Modify code so R is not computed unless needed
11 
16 #ifndef ANASAZI_TRACEMIN_BASE_HPP
17 #define ANASAZI_TRACEMIN_BASE_HPP
18 
19 #include "AnasaziBasicSort.hpp"
20 #include "AnasaziConfigDefs.hpp"
21 #include "AnasaziEigensolver.hpp"
27 #include "AnasaziSolverUtils.hpp"
29 #include "AnasaziTraceMinTypes.hpp"
30 #include "AnasaziTypes.hpp"
31 
33 #include "Teuchos_ScalarTraits.hpp"
36 #include "Teuchos_TimeMonitor.hpp"
37 
38 using Teuchos::RCP;
39 using Teuchos::rcp;
40 
41 namespace Anasazi {
52 namespace Experimental {
53 
55 
56 
61  template <class ScalarType, class MV>
64  int curDim;
93  bool isOrtho;
95  int NEV;
97  ScalarType largestSafeShift;
100  TraceMinBaseState() : curDim(0), V(Teuchos::null), KV(Teuchos::null), MopV(Teuchos::null),
101  X(Teuchos::null), KX(Teuchos::null), MX(Teuchos::null), R(Teuchos::null),
102  T(Teuchos::null), KK(Teuchos::null), RV(Teuchos::null), isOrtho(false),
103  NEV(0), largestSafeShift(Teuchos::ScalarTraits<ScalarType>::zero()),
104  ritzShifts(Teuchos::null) {}
105  };
106 
108 
110 
111 
124  class TraceMinBaseInitFailure : public AnasaziError {public:
125  TraceMinBaseInitFailure(const std::string& what_arg) : AnasaziError(what_arg)
126  {}};
127 
131  class TraceMinBaseOrthoFailure : public AnasaziError {public:
132  TraceMinBaseOrthoFailure(const std::string& what_arg) : AnasaziError(what_arg)
133  {}};
134 
136 
149  template <class ScalarType, class MV, class OP>
150  class TraceMinBase : public Eigensolver<ScalarType,MV,OP> {
151  public:
153 
154 
184  const RCP<OutputManager<ScalarType> > &printer,
185  const RCP<StatusTest<ScalarType,MV,OP> > &tester,
186  const RCP<MatOrthoManager<ScalarType,MV,OP> > &ortho,
187  Teuchos::ParameterList &params
188  );
189 
191  virtual ~TraceMinBase();
193 
194 
196 
197 
221  void iterate();
222 
223  void harmonicIterate();
224 
248 
249  void harmonicInitialize(TraceMinBaseState<ScalarType,MV> newstate);
250 
254  void initialize();
255 
271  bool isInitialized() const;
272 
282 
284 
285 
287 
288 
290  int getNumIters() const;
291 
293  void resetNumIters();
294 
303 
309  std::vector<Value<ScalarType> > getRitzValues();
310 
311 
320  std::vector<int> getRitzIndex();
321 
322 
328  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getResNorms();
329 
330 
336  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getRes2Norms();
337 
338 
346  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getRitzRes2Norms();
347 
353  int getCurSubspaceDim() const;
354 
356  int getMaxSubspaceDim() const;
357 
359 
360 
362 
363 
366 
369 
372 
382  void setBlockSize(int blockSize);
383 
385  int getBlockSize() const;
386 
404  void setAuxVecs(const Teuchos::Array<RCP<const MV> > &auxvecs);
405 
408 
410 
412 
413 
420  void setSize(int blockSize, int numBlocks);
421 
423 
424 
426  void currentStatus(std::ostream &os);
427 
429 
430  protected:
431  //
432  // Convenience typedefs
433  //
438  typedef typename SCT::magnitudeType MagnitudeType;
439  typedef TraceMinRitzOp<ScalarType,MV,OP> tracemin_ritz_op_type;
440  typedef SaddleContainer<ScalarType,MV> saddle_container_type;
441  typedef SaddleOperator<ScalarType,MV,tracemin_ritz_op_type> saddle_op_type;
442  const MagnitudeType ONE;
443  const MagnitudeType ZERO;
444  const MagnitudeType NANVAL;
445  //
446  // Classes inputed through constructor that define the eigenproblem to be solved.
447  //
448  const RCP<Eigenproblem<ScalarType,MV,OP> > problem_;
449  const RCP<SortManager<MagnitudeType> > sm_;
450  const RCP<OutputManager<ScalarType> > om_;
452  const RCP<MatOrthoManager<ScalarType,MV,OP> > orthman_;
453  //
454  // Information obtained from the eigenproblem
455  //
456  RCP<const OP> Op_;
457  RCP<const OP> MOp_;
458  RCP<const OP> Prec_;
459  bool hasM_;
460  //
461  // Internal timers
462  // TODO: Fix the timers
463  //
464  RCP<Teuchos::Time> timerOp_, timerMOp_, timerSaddle_, timerSortEval_, timerDS_,
465  timerLocal_, timerCompRes_, timerOrtho_, timerInit_;
466  //
467  // Internal structs
468  // TODO: Fix the checks
469  //
470  struct CheckList {
471  bool checkV, checkX, checkMX,
472  checkKX, checkQ, checkKK;
473  CheckList() : checkV(false),checkX(false),
474  checkMX(false),checkKX(false),
475  checkQ(false),checkKK(false) {};
476  };
477  //
478  // Internal methods
479  //
480  std::string accuracyCheck(const CheckList &chk, const std::string &where) const;
481 
482  //
483  // Counters
484  //
485  int count_ApplyOp_, count_ApplyM_;
486 
487  //
488  // Algorithmic parameters.
489  //
490  // blockSize_ is the solver block size; it controls the number of vectors added to the basis on each iteration.
491  int blockSize_;
492  // numBlocks_ is the size of the allocated space for the basis, in blocks.
493  int numBlocks_;
494 
495  //
496  // Current solver state
497  //
498  // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
499  // is capable of running; _initialize is controlled by the initialize() member method
500  // For the implications of the state of initialized_, please see documentation for initialize()
501  bool initialized_;
502  //
503  // curDim_ reflects how much of the current basis is valid
504  // NOTE: 0 <= curDim_ <= blockSize_*numBlocks_
505  // this also tells us how many of the values in theta_ are valid Ritz values
506  int curDim_;
507  //
508  // State Multivecs
509  // V is the current basis
510  // X is the current set of eigenvectors
511  // R is the residual
512  RCP<MV> X_, KX_, MX_, KV_, MV_, R_, V_;
513  //
514  // Projected matrices
515  //
517  //
518  // auxiliary vectors
519  Teuchos::Array<RCP<const MV> > auxVecs_, MauxVecs_;
520  int numAuxVecs_;
521  //
522  // Number of iterations that have been performed.
523  int iter_;
524  //
525  // Current eigenvalues, residual norms
526  std::vector<MagnitudeType> theta_, Rnorms_, R2norms_;
527  //
528  // are the residual norms current with the residual?
529  bool Rnorms_current_, R2norms_current_;
530 
531  // TraceMin specific variables
532  RCP<tracemin_ritz_op_type> ritzOp_; // Operator which incorporates Ritz shifts
533  enum SaddleSolType saddleSolType_; // Type of saddle point solver being used
534  bool previouslyLeveled_; // True if the trace already leveled off
535  MagnitudeType previousTrace_; // Trace of X'KX at the previous iteration
536  bool posSafeToShift_, negSafeToShift_; // Whether there are multiple clusters
537  MagnitudeType largestSafeShift_; // The largest shift considered to be safe - generally the biggest converged eigenvalue
538  int NEV_; // Number of eigenvalues we seek - used in computation of trace
539  std::vector<ScalarType> ritzShifts_; // Current Ritz shifts
540 
541  // This is only used if we're using the Schur complement solver
542  RCP<MV> Z_;
543 
544  // User provided TraceMin parameters
545  enum WhenToShiftType whenToShift_; // What triggers a Ritz shift
546  enum HowToShiftType howToShift_; // Strategy for choosing the Ritz shifts
547  bool useMultipleShifts_; // Whether to use one Ritz shift or many
548  bool considerClusters_; // Don't shift if there is only one cluster
549  bool projectAllVecs_; // Use all vectors in projected minres, or just 1
550  bool projectLockedVecs_; // Project locked vectors too in minres; does nothing if projectAllVecs = false
551  bool computeAllRes_; // Compute all residuals, or just blocksize ones - helps make Ritz shifts safer
552  bool useRHSR_; // Use -R as the RHS of projected minres rather than AX
553  bool useHarmonic_;
554  MagnitudeType traceThresh_;
555  MagnitudeType alpha_;
556 
557  // Variables that are only used if we're shifting when the residual gets small
558  // TODO: These are currently unused
559  std::string shiftNorm_; // Measure 2-norm or M-norm of residual
560  MagnitudeType shiftThresh_; // How small must the residual be?
561  bool useRelShiftThresh_; // Are we scaling the shift threshold by the eigenvalues?
562 
563  // TraceMin specific functions
564  // Returns the trace of KK = X'KX
565  ScalarType getTrace() const;
566  // Returns true if the change in trace is very small (or was very small at one point)
567  // TODO: Check whether I want to return true if the trace WAS small
568  bool traceLeveled();
569  // Get the residuals of each cluster of eigenvalues
570  // TODO: Figure out whether I want to use these for all eigenvalues or just the first
571  std::vector<ScalarType> getClusterResids();
572  // Computes the Ritz shifts, which is not a trivial task
573  // TODO: Make it safer for indefinite matrices maybe?
574  void computeRitzShifts(const std::vector<ScalarType>& clusterResids);
575  // Compute the tolerance for the inner solves
576  // TODO: Come back to this and make sure it does what I want it to
577  std::vector<ScalarType> computeTol();
578  // Solves a saddle point problem
579  void solveSaddlePointProblem(RCP<MV> Delta);
580  // Solves a saddle point problem by using unpreconditioned projected minres
581  void solveSaddleProj(RCP<MV> Delta) const;
582  // Solves a saddle point problem by using preconditioned projected...Krylov...something
583  // TODO: Fix this. Replace minres with gmres or something.
584  void solveSaddleProjPrec(RCP<MV> Delta) const;
585  // Solves a saddle point problem by explicitly forming the inexact Schur complement
586  void solveSaddleSchur (RCP<MV> Delta) const;
587  // Solves a saddle point problem with a block diagonal preconditioner
588  void solveSaddleBDPrec (RCP<MV> Delta) const;
589  // Solves a saddle point problem with a Hermitian/non-Hermitian splitting preconditioner
590  void solveSaddleHSSPrec (RCP<MV> Delta) const;
591  // Computes KK = X'KX
592  void computeKK();
593  // Computes the eigenpairs of KK
594  void computeRitzPairs();
595  // Computes X = V evecs
596  void computeX();
597  // Updates KX and MX without a matvec by MAGIC (or by multiplying KV and MV by evecs)
598  void updateKXMX();
599  // Updates the residual
600  void updateResidual();
601  // Adds Delta to the basis
602  // TraceMin and TraceMinDavidson each do this differently, which is why this is pure virtual
603  virtual void addToBasis(const RCP<const MV> Delta) =0;
604  virtual void harmonicAddToBasis(const RCP<const MV> Delta) =0;
605  };
606 
609  //
610  // Implementations
611  //
614 
615 
617  // Constructor
618  // TODO: Allow the users to supply their own linear solver
619  // TODO: Add additional checking for logic errors (like trying to use gmres with multiple shifts)
620  template <class ScalarType, class MV, class OP>
622  const RCP<Eigenproblem<ScalarType,MV,OP> > &problem,
624  const RCP<OutputManager<ScalarType> > &printer,
625  const RCP<StatusTest<ScalarType,MV,OP> > &tester,
626  const RCP<MatOrthoManager<ScalarType,MV,OP> > &ortho,
627  Teuchos::ParameterList &params
628  ) :
629  ONE(Teuchos::ScalarTraits<MagnitudeType>::one()),
630  ZERO(Teuchos::ScalarTraits<MagnitudeType>::zero()),
631  NANVAL(Teuchos::ScalarTraits<MagnitudeType>::nan()),
632  // problem, tools
633  problem_(problem),
634  sm_(sorter),
635  om_(printer),
636  tester_(tester),
637  orthman_(ortho),
638  // timers, counters
639 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
640  timerOp_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Operation Op*x")),
641  timerMOp_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Operation M*x")),
642  timerSaddle_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Solving saddle point problem")),
643  timerSortEval_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Sorting eigenvalues")),
644  timerDS_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Direct solve")),
645  timerLocal_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Local update")),
646  timerCompRes_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Computing residuals")),
647  timerOrtho_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Orthogonalization")),
648  timerInit_(Teuchos::TimeMonitor::getNewTimer("Anasazi: TraceMinBase::Initialization")),
649 #endif
650  count_ApplyOp_(0),
651  count_ApplyM_(0),
652  // internal data
653  blockSize_(0),
654  numBlocks_(0),
655  initialized_(false),
656  curDim_(0),
657  auxVecs_( Teuchos::Array<RCP<const MV> >(0) ),
658  MauxVecs_( Teuchos::Array<RCP<const MV> >(0) ),
659  numAuxVecs_(0),
660  iter_(0),
661  Rnorms_current_(false),
662  R2norms_current_(false),
663  // TraceMin variables
664  previouslyLeveled_(false),
665  previousTrace_(ZERO),
666  posSafeToShift_(false),
667  negSafeToShift_(false),
668  largestSafeShift_(ZERO),
669  Z_(Teuchos::null)
670  {
671  TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null,std::invalid_argument,
672  "Anasazi::TraceMinBase::constructor: user passed null problem pointer.");
673  TEUCHOS_TEST_FOR_EXCEPTION(sm_ == Teuchos::null,std::invalid_argument,
674  "Anasazi::TraceMinBase::constructor: user passed null sort manager pointer.");
675  TEUCHOS_TEST_FOR_EXCEPTION(om_ == Teuchos::null,std::invalid_argument,
676  "Anasazi::TraceMinBase::constructor: user passed null output manager pointer.");
677  TEUCHOS_TEST_FOR_EXCEPTION(tester_ == Teuchos::null,std::invalid_argument,
678  "Anasazi::TraceMinBase::constructor: user passed null status test pointer.");
679  TEUCHOS_TEST_FOR_EXCEPTION(orthman_ == Teuchos::null,std::invalid_argument,
680  "Anasazi::TraceMinBase::constructor: user passed null orthogonalization manager pointer.");
681  TEUCHOS_TEST_FOR_EXCEPTION(problem_->isHermitian() == false, std::invalid_argument,
682  "Anasazi::TraceMinBase::constructor: problem is not hermitian.");
683 
684  // get the problem operators
685  Op_ = problem_->getOperator();
686  MOp_ = problem_->getM();
687  Prec_ = problem_->getPrec();
688  hasM_ = (MOp_ != Teuchos::null);
689 
690  // Set the saddle point solver parameters
691  saddleSolType_ = params.get("Saddle Solver Type", PROJECTED_KRYLOV_SOLVER);
692  TEUCHOS_TEST_FOR_EXCEPTION(saddleSolType_ != PROJECTED_KRYLOV_SOLVER && saddleSolType_ != SCHUR_COMPLEMENT_SOLVER && saddleSolType_ != BD_PREC_MINRES && saddleSolType_ != HSS_PREC_GMRES, std::invalid_argument,
693  "Anasazi::TraceMin::constructor: Invalid value for \"Saddle Solver Type\"; valid options are PROJECTED_KRYLOV_SOLVER, SCHUR_COMPLEMENT_SOLVER, and BD_PREC_MINRES.");
694 
695  // Set the Ritz shift parameters
696  whenToShift_ = params.get("When To Shift", ALWAYS_SHIFT);
697  TEUCHOS_TEST_FOR_EXCEPTION(whenToShift_ != NEVER_SHIFT && whenToShift_ != SHIFT_WHEN_TRACE_LEVELS && whenToShift_ != SHIFT_WHEN_RESID_SMALL && whenToShift_ != ALWAYS_SHIFT, std::invalid_argument,
698  "Anasazi::TraceMin::constructor: Invalid value for \"When To Shift\"; valid options are \"NEVER_SHIFT\", \"SHIFT_WHEN_TRACE_LEVELS\", \"SHIFT_WHEN_RESID_SMALL\", and \"ALWAYS_SHIFT\".");
699 
700  traceThresh_ = params.get("Trace Threshold", 2e-2);
701  TEUCHOS_TEST_FOR_EXCEPTION(traceThresh_ < 0, std::invalid_argument,
702  "Anasazi::TraceMin::constructor: Invalid value for \"Trace Threshold\"; Must be positive.");
703 
704  howToShift_ = params.get("How To Choose Shift", ADJUSTED_RITZ_SHIFT);
705  TEUCHOS_TEST_FOR_EXCEPTION(howToShift_ != LARGEST_CONVERGED_SHIFT && howToShift_ != ADJUSTED_RITZ_SHIFT && howToShift_ != RITZ_VALUES_SHIFT && howToShift_ != EXPERIMENTAL_SHIFT, std::invalid_argument,
706  "Anasazi::TraceMin::constructor: Invalid value for \"How To Choose Shift\"; valid options are \"LARGEST_CONVERGED_SHIFT\", \"ADJUSTED_RITZ_SHIFT\", \"RITZ_VALUES_SHIFT\".");
707 
708  useMultipleShifts_ = params.get("Use Multiple Shifts", true);
709 
710  shiftThresh_ = params.get("Shift Threshold", 1e-2);
711  useRelShiftThresh_ = params.get("Relative Shift Threshold", true);
712  shiftNorm_ = params.get("Shift Norm", "2");
713  TEUCHOS_TEST_FOR_EXCEPTION(shiftNorm_ != "2" && shiftNorm_ != "M", std::invalid_argument,
714  "Anasazi::TraceMin::constructor: Invalid value for \"Shift Norm\"; valid options are \"2\", \"M\".");
715 
716  considerClusters_ = params.get("Consider Clusters", true);
717 
718  projectAllVecs_ = params.get("Project All Vectors", true);
719  projectLockedVecs_ = params.get("Project Locked Vectors", true);
720  useRHSR_ = params.get("Use Residual as RHS", true);
721  useHarmonic_ = params.get("Use Harmonic Ritz Values", false);
722  computeAllRes_ = params.get("Compute All Residuals", true);
723 
724  // set the block size and allocate data
725  int bs = params.get("Block Size", problem_->getNEV());
726  int nb = params.get("Num Blocks", 1);
727  setSize(bs,nb);
728 
729  NEV_ = problem_->getNEV();
730 
731  // Create the Ritz shift operator
732  ritzOp_ = rcp (new tracemin_ritz_op_type (Op_, MOp_, Prec_));
733 
734  // Set the maximum number of inner iterations
735  const int innerMaxIts = params.get ("Maximum Krylov Iterations", 200);
736  ritzOp_->setMaxIts (innerMaxIts);
737 
738  alpha_ = params.get ("HSS: alpha", ONE);
739  }
740 
741 
743  // Destructor
744  template <class ScalarType, class MV, class OP>
746 
747 
749  // Set the block size
750  // This simply calls setSize(), modifying the block size while retaining the number of blocks.
751  template <class ScalarType, class MV, class OP>
753  {
754  TEUCHOS_TEST_FOR_EXCEPTION(blockSize < 1, std::invalid_argument, "Anasazi::TraceMinBase::setSize(blocksize,numblocks): blocksize must be strictly positive.");
755  setSize(blockSize,numBlocks_);
756  }
757 
758 
760  // Return the current auxiliary vectors
761  template <class ScalarType, class MV, class OP>
763  return auxVecs_;
764  }
765 
766 
768  // return the current block size
769  template <class ScalarType, class MV, class OP>
771  return(blockSize_);
772  }
773 
774 
776  // return eigenproblem
777  template <class ScalarType, class MV, class OP>
779  return(*problem_);
780  }
781 
782 
784  // return max subspace dim
785  template <class ScalarType, class MV, class OP>
787  return blockSize_*numBlocks_;
788  }
789 
791  // return current subspace dim
792  template <class ScalarType, class MV, class OP>
794  if (!initialized_) return 0;
795  return curDim_;
796  }
797 
798 
800  // return ritz residual 2-norms
801  template <class ScalarType, class MV, class OP>
802  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
804  return getRes2Norms();
805  }
806 
807 
809  // return ritz index
810  template <class ScalarType, class MV, class OP>
812  std::vector<int> ret(curDim_,0);
813  return ret;
814  }
815 
816 
818  // return ritz values
819  template <class ScalarType, class MV, class OP>
820  std::vector<Value<ScalarType> > TraceMinBase<ScalarType,MV,OP>::getRitzValues() {
821  std::vector<Value<ScalarType> > ret(curDim_);
822  for (int i=0; i<curDim_; ++i) {
823  ret[i].realpart = theta_[i];
824  ret[i].imagpart = ZERO;
825  }
826  return ret;
827  }
828 
829 
831  // return pointer to ritz vectors
832  template <class ScalarType, class MV, class OP>
834  return X_;
835  }
836 
837 
839  // reset number of iterations
840  template <class ScalarType, class MV, class OP>
842  iter_=0;
843  }
844 
845 
847  // return number of iterations
848  template <class ScalarType, class MV, class OP>
850  return(iter_);
851  }
852 
853 
855  // return state pointers
856  template <class ScalarType, class MV, class OP>
859  state.curDim = curDim_;
860  state.V = V_;
861  state.X = X_;
862  if (Op_ != Teuchos::null) {
863  state.KV = KV_;
864  state.KX = KX_;
865  }
866  else {
867  state.KV = Teuchos::null;
868  state.KX = Teuchos::null;
869  }
870  if (hasM_) {
871  state.MopV = MV_;
872  state.MX = MX_;
873  }
874  else {
875  state.MopV = Teuchos::null;
876  state.MX = Teuchos::null;
877  }
878  state.R = R_;
879  state.KK = KK_;
880  state.RV = ritzVecs_;
881  if (curDim_ > 0) {
882  state.T = rcp(new std::vector<MagnitudeType>(theta_.begin(),theta_.begin()+curDim_) );
883  } else {
884  state.T = rcp(new std::vector<MagnitudeType>(0));
885  }
886  state.ritzShifts = rcp(new std::vector<MagnitudeType>(ritzShifts_.begin(),ritzShifts_.begin()+blockSize_) );
887  state.isOrtho = true;
888  state.largestSafeShift = largestSafeShift_;
889 
890  return state;
891  }
892 
893 
895  // Perform TraceMinBase iterations until the StatusTest tells us to stop.
896  template <class ScalarType, class MV, class OP>
898  {
899  if(useHarmonic_)
900  {
901  harmonicIterate();
902  return;
903  }
904 
905  //
906  // Initialize solver state
907  if (initialized_ == false) {
908  initialize();
909  }
910 
911  // as a data member, this would be redundant and require synchronization with
912  // blockSize_ and numBlocks_; we'll use a constant here.
913  const int searchDim = blockSize_*numBlocks_;
914 
915  // Update obtained from solving saddle point problem
916  RCP<MV> Delta = MVT::Clone(*X_,blockSize_);
917 
919  // iterate until the status test tells us to stop.
920  // also break if our basis is full
921  while (tester_->checkStatus(this) != Passed && (numBlocks_ == 1 || curDim_ < searchDim)) {
922 
923  // Print information on current iteration
924  if (om_->isVerbosity(Debug)) {
925  currentStatus( om_->stream(Debug) );
926  }
927  else if (om_->isVerbosity(IterationDetails)) {
928  currentStatus( om_->stream(IterationDetails) );
929  }
930 
931  ++iter_;
932 
933  // Solve the saddle point problem
934  solveSaddlePointProblem(Delta);
935 
936  // Insert Delta at the end of V
937  addToBasis(Delta);
938 
939  if (om_->isVerbosity( Debug ) ) {
940  // Check almost everything here
941  CheckList chk;
942  chk.checkV = true;
943  om_->print( Debug, accuracyCheck(chk, ": after addToBasis(Delta)") );
944  }
945 
946  // Compute KK := V'KV
947  computeKK();
948 
949  if (om_->isVerbosity( Debug ) ) {
950  // Check almost everything here
951  CheckList chk;
952  chk.checkKK = true;
953  om_->print( Debug, accuracyCheck(chk, ": after computeKK()") );
954  }
955 
956  // Compute the Ritz pairs
957  computeRitzPairs();
958 
959  // Compute X := V RitzVecs
960  computeX();
961 
962  if (om_->isVerbosity( Debug ) ) {
963  // Check almost everything here
964  CheckList chk;
965  chk.checkX = true;
966  om_->print( Debug, accuracyCheck(chk, ": after computeX()") );
967  }
968 
969  // Compute KX := KV RitzVecs and MX := MV RitzVecs (if necessary)
970  updateKXMX();
971 
972  if (om_->isVerbosity( Debug ) ) {
973  // Check almost everything here
974  CheckList chk;
975  chk.checkKX = true;
976  chk.checkMX = true;
977  om_->print( Debug, accuracyCheck(chk, ": after updateKXMX()") );
978  }
979 
980  // Update the residual vectors
981  updateResidual();
982  } // end while (statusTest == false)
983 
984  } // end of iterate()
985 
986 
987 
989  // Perform TraceMinBase iterations until the StatusTest tells us to stop.
990  template <class ScalarType, class MV, class OP>
992  {
993  //
994  // Initialize solver state
995  if (initialized_ == false) {
996  initialize();
997  }
998 
999  // as a data member, this would be redundant and require synchronization with
1000  // blockSize_ and numBlocks_; we'll use a constant here.
1001  const int searchDim = blockSize_*numBlocks_;
1002 
1003  // Update obtained from solving saddle point problem
1004  RCP<MV> Delta = MVT::Clone(*X_,blockSize_);
1005 
1007  // iterate until the status test tells us to stop.
1008  // also break if our basis is full
1009  while (tester_->checkStatus(this) != Passed && (numBlocks_ == 1 || curDim_ < searchDim)) {
1010 
1011  // Print information on current iteration
1012  if (om_->isVerbosity(Debug)) {
1013  currentStatus( om_->stream(Debug) );
1014  }
1015  else if (om_->isVerbosity(IterationDetails)) {
1016  currentStatus( om_->stream(IterationDetails) );
1017  }
1018 
1019  ++iter_;
1020 
1021  // Solve the saddle point problem
1022  solveSaddlePointProblem(Delta);
1023 
1024  // Insert Delta at the end of V
1025  harmonicAddToBasis(Delta);
1026 
1027  if (om_->isVerbosity( Debug ) ) {
1028  // Check almost everything here
1029  CheckList chk;
1030  chk.checkV = true;
1031  om_->print( Debug, accuracyCheck(chk, ": after addToBasis(Delta)") );
1032  }
1033 
1034  // Compute KK := V'KV
1035  computeKK();
1036 
1037  if (om_->isVerbosity( Debug ) ) {
1038  // Check almost everything here
1039  CheckList chk;
1040  chk.checkKK = true;
1041  om_->print( Debug, accuracyCheck(chk, ": after computeKK()") );
1042  }
1043 
1044  // Compute the Ritz pairs
1045  computeRitzPairs();
1046 
1047  // Compute X := V RitzVecs
1048  computeX();
1049 
1050  // Get norm of each vector in X
1051  int nvecs;
1052  if(computeAllRes_)
1053  nvecs = curDim_;
1054  else
1055  nvecs = blockSize_;
1056  std::vector<int> dimind(nvecs);
1057  for(int i=0; i<nvecs; i++)
1058  dimind[i] = i;
1059  RCP<MV> lclX = MVT::CloneViewNonConst(*X_,dimind);
1060  std::vector<ScalarType> normvec(nvecs);
1061  orthman_->normMat(*lclX,normvec);
1062 
1063  // Scale X
1064  for(int i=0; i<nvecs; i++)
1065  normvec[i] = ONE/normvec[i];
1066  MVT::MvScale(*lclX,normvec);
1067 
1068  // Scale eigenvalues
1069  for(int i=0; i<nvecs; i++)
1070  {
1071  theta_[i] = theta_[i] * normvec[i] * normvec[i];
1072  }
1073 
1074  if (om_->isVerbosity( Debug ) ) {
1075  // Check almost everything here
1076  CheckList chk;
1077  chk.checkX = true;
1078  om_->print( Debug, accuracyCheck(chk, ": after computeX()") );
1079  }
1080 
1081  // Compute KX := KV RitzVecs and MX := MV RitzVecs (if necessary)
1082  updateKXMX();
1083 
1084  // Scale KX and MX
1085  if(Op_ != Teuchos::null)
1086  {
1087  RCP<MV> lclKX = MVT::CloneViewNonConst(*KX_,dimind);
1088  MVT::MvScale(*lclKX,normvec);
1089  }
1090  if(hasM_)
1091  {
1092  RCP<MV> lclMX = MVT::CloneViewNonConst(*MX_,dimind);
1093  MVT::MvScale(*lclMX,normvec);
1094  }
1095 
1096  if (om_->isVerbosity( Debug ) ) {
1097  // Check almost everything here
1098  CheckList chk;
1099  chk.checkKX = true;
1100  chk.checkMX = true;
1101  om_->print( Debug, accuracyCheck(chk, ": after updateKXMX()") );
1102  }
1103 
1104  // Update the residual vectors
1105  updateResidual();
1106  } // end while (statusTest == false)
1107 
1108  } // end of harmonicIterate()
1109 
1110 
1112  // initialize the solver with default state
1113  template <class ScalarType, class MV, class OP>
1115  {
1117  initialize(empty);
1118  }
1119 
1120 
1122  /* Initialize the state of the solver
1123  *
1124  * POST-CONDITIONS:
1125  *
1126  * V_ is orthonormal, orthogonal to auxVecs_, for first curDim_ vectors
1127  * theta_ contains Ritz w.r.t. V_(1:curDim_)
1128  * X is Ritz vectors w.r.t. V_(1:curDim_)
1129  * KX = Op*X
1130  * MX = M*X if hasM_
1131  * R = KX - MX*diag(theta_)
1132  *
1133  */
1134  template <class ScalarType, class MV, class OP>
1136  {
1137  // TODO: Check for bad input
1138  // NOTE: memory has been allocated by setBlockSize(). Use setBlock below; do not Clone
1139  // NOTE: Overall time spent in this routine is counted to timerInit_; portions will also be counted towards other primitives
1140 
1141 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1142  Teuchos::TimeMonitor inittimer( *timerInit_ );
1143 #endif
1144 
1145  previouslyLeveled_ = false;
1146 
1147  if(useHarmonic_)
1148  {
1149  harmonicInitialize(newstate);
1150  return;
1151  }
1152 
1153  std::vector<int> bsind(blockSize_);
1154  for (int i=0; i<blockSize_; ++i) bsind[i] = i;
1155 
1156  // in TraceMinBase, V is primary
1157  // the order of dependence follows like so.
1158  // --init-> V,KK
1159  // --ritz analysis-> theta,X
1160  // --op apply-> KX,MX
1161  // --compute-> R
1162  //
1163  // if the user specifies all data for a level, we will accept it.
1164  // otherwise, we will generate the whole level, and all subsequent levels.
1165  //
1166  // the data members are ordered based on dependence, and the levels are
1167  // partitioned according to the amount of work required to produce the
1168  // items in a level.
1169  //
1170  // inconsistent multivectors widths and lengths will not be tolerated, and
1171  // will be treated with exceptions.
1172  //
1173  // for multivector pointers in newstate which point directly (as opposed to indirectly, via a view) to
1174  // multivectors in the solver, the copy will not be affected.
1175 
1176  // set up V and KK: get them from newstate if user specified them
1177  // otherwise, set them manually
1178  RCP<MV> lclV, lclKV, lclMV;
1180 
1181  // If V was supplied by the user...
1182  if (newstate.V != Teuchos::null) {
1183  om_->stream(Debug) << "Copying V from the user\n";
1184 
1185  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.V) != MVT::GetGlobalLength(*V_), std::invalid_argument,
1186  "Anasazi::TraceMinBase::initialize(newstate): Vector length of V not correct." );
1187  TEUCHOS_TEST_FOR_EXCEPTION( newstate.curDim < blockSize_, std::invalid_argument,
1188  "Anasazi::TraceMinBase::initialize(newstate): Rank of new state must be at least blockSize().");
1189  TEUCHOS_TEST_FOR_EXCEPTION( newstate.curDim > blockSize_*numBlocks_, std::invalid_argument,
1190  "Anasazi::TraceMinBase::initialize(newstate): Rank of new state must be less than getMaxSubspaceDim().");
1191 
1192  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.V) < newstate.curDim, std::invalid_argument,
1193  "Anasazi::TraceMinBase::initialize(newstate): Multivector for basis in new state must be as large as specified state rank.");
1194 
1195  curDim_ = newstate.curDim;
1196  // pick an integral amount
1197  curDim_ = (int)(curDim_ / blockSize_)*blockSize_;
1198 
1199  TEUCHOS_TEST_FOR_EXCEPTION( curDim_ != newstate.curDim, std::invalid_argument,
1200  "Anasazi::TraceMinBase::initialize(newstate): Rank of new state must be a multiple of getBlockSize().");
1201 
1202  // put data in V
1203  std::vector<int> nevind(curDim_);
1204  for (int i=0; i<curDim_; ++i) nevind[i] = i;
1205  if (newstate.V != V_) {
1206  if (curDim_ < MVT::GetNumberVecs(*newstate.V)) {
1207  newstate.V = MVT::CloneView(*newstate.V,nevind);
1208  }
1209  MVT::SetBlock(*newstate.V,nevind,*V_);
1210  }
1211  lclV = MVT::CloneViewNonConst(*V_,nevind);
1212  } // end if user specified V
1213  else
1214  {
1215  // get vectors from problem or generate something, projectAndNormalize
1216  RCP<const MV> ivec = problem_->getInitVec();
1217  TEUCHOS_TEST_FOR_EXCEPTION(ivec == Teuchos::null,std::invalid_argument,
1218  "Anasazi::TraceMinBase::initialize(newstate): Eigenproblem did not specify initial vectors to clone from.");
1219  // clear newstate so we won't use any data from it below
1220  newstate.X = Teuchos::null;
1221  newstate.MX = Teuchos::null;
1222  newstate.KX = Teuchos::null;
1223  newstate.R = Teuchos::null;
1224  newstate.T = Teuchos::null;
1225  newstate.RV = Teuchos::null;
1226  newstate.KK = Teuchos::null;
1227  newstate.KV = Teuchos::null;
1228  newstate.MopV = Teuchos::null;
1229  newstate.isOrtho = false;
1230 
1231  // If the user did not specify a current dimension, use that of the initial vectors they provided
1232  if(newstate.curDim > 0)
1233  curDim_ = newstate.curDim;
1234  else
1235  curDim_ = MVT::GetNumberVecs(*ivec);
1236 
1237  // pick the largest multiple of blockSize_
1238  curDim_ = (int)(curDim_ / blockSize_)*blockSize_;
1239  if (curDim_ > blockSize_*numBlocks_) {
1240  // user specified too many vectors... truncate
1241  // this produces a full subspace, but that is okay
1242  curDim_ = blockSize_*numBlocks_;
1243  }
1244  bool userand = false;
1245  if (curDim_ == 0) {
1246  // we need at least blockSize_ vectors
1247  // use a random multivec: ignore everything from InitVec
1248  userand = true;
1249  curDim_ = blockSize_;
1250  }
1251 
1252  std::vector<int> nevind(curDim_);
1253  for (int i=0; i<curDim_; ++i) nevind[i] = i;
1254 
1255  // get a pointer into V
1256  // lclV has curDim vectors
1257  //
1258  // get pointer to first curDim vectors in V_
1259  lclV = MVT::CloneViewNonConst(*V_,nevind);
1260  if (userand)
1261  {
1262  // generate random vector data
1263  MVT::MvRandom(*lclV);
1264  }
1265  else
1266  {
1267  if(newstate.curDim > 0)
1268  {
1269  if(MVT::GetNumberVecs(*newstate.V) > curDim_) {
1270  RCP<const MV> helperMV = MVT::CloneView(*newstate.V,nevind);
1271  MVT::SetBlock(*helperMV,nevind,*lclV);
1272  }
1273  // assign ivec to first part of lclV
1274  MVT::SetBlock(*newstate.V,nevind,*lclV);
1275  }
1276  else
1277  {
1278  if(MVT::GetNumberVecs(*ivec) > curDim_) {
1279  ivec = MVT::CloneView(*ivec,nevind);
1280  }
1281  // assign ivec to first part of lclV
1282  MVT::SetBlock(*ivec,nevind,*lclV);
1283  }
1284  }
1285  } // end if user did not specify V
1286 
1287  // A vector of relevant indices
1288  std::vector<int> dimind(curDim_);
1289  for (int i=0; i<curDim_; ++i) dimind[i] = i;
1290 
1291  // Compute MV if necessary
1292  if(hasM_ && newstate.MopV == Teuchos::null)
1293  {
1294  om_->stream(Debug) << "Computing MV\n";
1295 
1296 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1297  Teuchos::TimeMonitor lcltimer( *timerMOp_ );
1298 #endif
1299  count_ApplyM_+= curDim_;
1300 
1301  newstate.isOrtho = false;
1302  // Get a pointer to the relevant parts of MV
1303  lclMV = MVT::CloneViewNonConst(*MV_,dimind);
1304  OPT::Apply(*MOp_,*lclV,*lclMV);
1305  }
1306  // Copy MV if necessary
1307  else if(hasM_)
1308  {
1309  om_->stream(Debug) << "Copying MV\n";
1310 
1311  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.MopV) != MVT::GetGlobalLength(*MV_), std::invalid_argument,
1312  "Anasazi::TraceMinBase::initialize(newstate): Vector length of MV not correct." );
1313 
1314  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.MopV) < curDim_, std::invalid_argument,
1315  "Anasazi::TraceMinBase::initialize(newstate): Number of vectors in MV not correct.");
1316 
1317  if(newstate.MopV != MV_) {
1318  if(curDim_ < MVT::GetNumberVecs(*newstate.MopV)) {
1319  newstate.MopV = MVT::CloneView(*newstate.MopV,dimind);
1320  }
1321  MVT::SetBlock(*newstate.MopV,dimind,*MV_);
1322  }
1323  lclMV = MVT::CloneViewNonConst(*MV_,dimind);
1324  }
1325  // There is no M, so there is no MV
1326  else
1327  {
1328  om_->stream(Debug) << "There is no MV\n";
1329 
1330  lclMV = lclV;
1331  MV_ = V_;
1332  }
1333 
1334  // Project and normalize so that V'V = I and Q'V=0
1335  if(newstate.isOrtho == false)
1336  {
1337  om_->stream(Debug) << "Project and normalize\n";
1338 
1339 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1340  Teuchos::TimeMonitor lcltimer( *timerOrtho_ );
1341 #endif
1342 
1343  // These things are now invalid
1344  newstate.X = Teuchos::null;
1345  newstate.MX = Teuchos::null;
1346  newstate.KX = Teuchos::null;
1347  newstate.R = Teuchos::null;
1348  newstate.T = Teuchos::null;
1349  newstate.RV = Teuchos::null;
1350  newstate.KK = Teuchos::null;
1351  newstate.KV = Teuchos::null;
1352 
1353  int rank;
1354  if(auxVecs_.size() > 0)
1355  {
1356  rank = orthman_->projectAndNormalizeMat(*lclV, auxVecs_,
1357  Teuchos::tuple(RCP< Teuchos::SerialDenseMatrix< int, ScalarType > >(Teuchos::null)),
1358  Teuchos::null, lclMV, MauxVecs_);
1359  }
1360  else
1361  {
1362  rank = orthman_->normalizeMat(*lclV,Teuchos::null,lclMV);
1363  }
1364 
1366  "Anasazi::TraceMinBase::initialize(): Couldn't generate initial basis of full rank.");
1367  }
1368 
1369  // Compute KV
1370  if(Op_ != Teuchos::null && newstate.KV == Teuchos::null)
1371  {
1372  om_->stream(Debug) << "Computing KV\n";
1373 
1374 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1375  Teuchos::TimeMonitor lcltimer( *timerOp_ );
1376 #endif
1377  count_ApplyOp_+= curDim_;
1378 
1379  // These things are now invalid
1380  newstate.X = Teuchos::null;
1381  newstate.MX = Teuchos::null;
1382  newstate.KX = Teuchos::null;
1383  newstate.R = Teuchos::null;
1384  newstate.T = Teuchos::null;
1385  newstate.RV = Teuchos::null;
1386  newstate.KK = Teuchos::null;
1387  newstate.KV = Teuchos::null;
1388 
1389  lclKV = MVT::CloneViewNonConst(*KV_,dimind);
1390  OPT::Apply(*Op_,*lclV,*lclKV);
1391  }
1392  // Copy KV
1393  else if(Op_ != Teuchos::null)
1394  {
1395  om_->stream(Debug) << "Copying MV\n";
1396 
1397  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.KV) != MVT::GetGlobalLength(*KV_), std::invalid_argument,
1398  "Anasazi::TraceMinBase::initialize(newstate): Vector length of MV not correct." );
1399 
1400  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.KV) < curDim_, std::invalid_argument,
1401  "Anasazi::TraceMinBase::initialize(newstate): Number of vectors in KV not correct.");
1402 
1403  if (newstate.KV != KV_) {
1404  if (curDim_ < MVT::GetNumberVecs(*newstate.KV)) {
1405  newstate.KV = MVT::CloneView(*newstate.KV,dimind);
1406  }
1407  MVT::SetBlock(*newstate.KV,dimind,*KV_);
1408  }
1409  lclKV = MVT::CloneViewNonConst(*KV_,dimind);
1410  }
1411  else
1412  {
1413  om_->stream(Debug) << "There is no KV\n";
1414 
1415  lclKV = lclV;
1416  KV_ = V_;
1417  }
1418 
1419  // Compute KK
1420  if(newstate.KK == Teuchos::null)
1421  {
1422  om_->stream(Debug) << "Computing KK\n";
1423 
1424  // These things are now invalid
1425  newstate.X = Teuchos::null;
1426  newstate.MX = Teuchos::null;
1427  newstate.KX = Teuchos::null;
1428  newstate.R = Teuchos::null;
1429  newstate.T = Teuchos::null;
1430  newstate.RV = Teuchos::null;
1431 
1432  // Note: there used to be a bug here.
1433  // We can't just call computeKK because it only computes the new parts of KK
1434 
1435  // Get a pointer to the part of KK we're interested in
1436  lclKK = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,curDim_) );
1437 
1438  // KK := V'KV
1439  MVT::MvTransMv(ONE,*lclV,*lclKV,*lclKK);
1440  }
1441  // Copy KK
1442  else
1443  {
1444  om_->stream(Debug) << "Copying KK\n";
1445 
1446  // check size of KK
1447  TEUCHOS_TEST_FOR_EXCEPTION( newstate.KK->numRows() < curDim_ || newstate.KK->numCols() < curDim_, std::invalid_argument,
1448  "Anasazi::TraceMinBase::initialize(newstate): Projected matrix in new state must be as large as specified state rank.");
1449 
1450  // put data into KK_
1451  lclKK = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,curDim_) );
1452  if (newstate.KK != KK_) {
1453  if (newstate.KK->numRows() > curDim_ || newstate.KK->numCols() > curDim_) {
1454  newstate.KK = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*newstate.KK,curDim_,curDim_) );
1455  }
1456  lclKK->assign(*newstate.KK);
1457  }
1458  }
1459 
1460  // Compute Ritz pairs
1461  if(newstate.T == Teuchos::null || newstate.RV == Teuchos::null)
1462  {
1463  om_->stream(Debug) << "Computing Ritz pairs\n";
1464 
1465  // These things are now invalid
1466  newstate.X = Teuchos::null;
1467  newstate.MX = Teuchos::null;
1468  newstate.KX = Teuchos::null;
1469  newstate.R = Teuchos::null;
1470  newstate.T = Teuchos::null;
1471  newstate.RV = Teuchos::null;
1472 
1473  computeRitzPairs();
1474  }
1475  // Copy Ritz pairs
1476  else
1477  {
1478  om_->stream(Debug) << "Copying Ritz pairs\n";
1479 
1480  TEUCHOS_TEST_FOR_EXCEPTION((signed int)(newstate.T->size()) != curDim_,
1481  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of T must be consistent with dimension of V.");
1482 
1483  TEUCHOS_TEST_FOR_EXCEPTION( newstate.RV->numRows() < curDim_ || newstate.RV->numCols() < curDim_, std::invalid_argument,
1484  "Anasazi::TraceMinBase::initialize(newstate): Ritz vectors in new state must be as large as specified state rank.");
1485 
1486  std::copy(newstate.T->begin(),newstate.T->end(),theta_.begin());
1487 
1488  lclRV = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,curDim_) );
1489  if (newstate.RV != ritzVecs_) {
1490  if (newstate.RV->numRows() > curDim_ || newstate.RV->numCols() > curDim_) {
1491  newstate.RV = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*newstate.RV,curDim_,curDim_) );
1492  }
1493  lclRV->assign(*newstate.RV);
1494  }
1495  }
1496 
1497  // Compute X
1498  if(newstate.X == Teuchos::null)
1499  {
1500  om_->stream(Debug) << "Computing X\n";
1501 
1502  // These things are now invalid
1503  newstate.MX = Teuchos::null;
1504  newstate.KX = Teuchos::null;
1505  newstate.R = Teuchos::null;
1506 
1507  computeX();
1508  }
1509  // Copy X
1510  else
1511  {
1512  om_->stream(Debug) << "Copying X\n";
1513 
1514  if(computeAllRes_ == false)
1515  {
1516  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.X) != blockSize_ || MVT::GetGlobalLength(*newstate.X) != MVT::GetGlobalLength(*X_),
1517  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of X must be consistent with block size and length of V.");
1518 
1519  if(newstate.X != X_) {
1520  MVT::SetBlock(*newstate.X,bsind,*X_);
1521  }
1522  }
1523  else
1524  {
1525  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.X) != curDim_ || MVT::GetGlobalLength(*newstate.X) != MVT::GetGlobalLength(*X_),
1526  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of X must be consistent with current dimension and length of V.");
1527 
1528  if(newstate.X != X_) {
1529  MVT::SetBlock(*newstate.X,dimind,*X_);
1530  }
1531  }
1532  }
1533 
1534  // Compute KX and MX if necessary
1535  // TODO: These technically should be separate; it won't matter much in terms of running time though
1536  if((Op_ != Teuchos::null && newstate.KX == Teuchos::null) || (hasM_ && newstate.MX == Teuchos::null))
1537  {
1538  om_->stream(Debug) << "Computing KX and MX\n";
1539 
1540  // These things are now invalid
1541  newstate.R = Teuchos::null;
1542 
1543  updateKXMX();
1544  }
1545  // Copy KX and MX if necessary
1546  else
1547  {
1548  om_->stream(Debug) << "Copying KX and MX\n";
1549  if(Op_ != Teuchos::null)
1550  {
1551  if(computeAllRes_ == false)
1552  {
1553  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.KX) != blockSize_ || MVT::GetGlobalLength(*newstate.KX) != MVT::GetGlobalLength(*KX_),
1554  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of KX must be consistent with block size and length of V.");
1555 
1556  if(newstate.KX != KX_) {
1557  MVT::SetBlock(*newstate.KX,bsind,*KX_);
1558  }
1559  }
1560  else
1561  {
1562  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.KX) != curDim_ || MVT::GetGlobalLength(*newstate.KX) != MVT::GetGlobalLength(*KX_),
1563  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of KX must be consistent with current dimension and length of V.");
1564 
1565  if (newstate.KX != KX_) {
1566  MVT::SetBlock(*newstate.KX,dimind,*KX_);
1567  }
1568  }
1569  }
1570 
1571  if(hasM_)
1572  {
1573  if(computeAllRes_ == false)
1574  {
1575  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.MX) != blockSize_ || MVT::GetGlobalLength(*newstate.MX) != MVT::GetGlobalLength(*MX_),
1576  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of MX must be consistent with block size and length of V.");
1577 
1578  if (newstate.MX != MX_) {
1579  MVT::SetBlock(*newstate.MX,bsind,*MX_);
1580  }
1581  }
1582  else
1583  {
1584  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.MX) != curDim_ || MVT::GetGlobalLength(*newstate.MX) != MVT::GetGlobalLength(*MX_),
1585  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of MX must be consistent with current dimension and length of V.");
1586 
1587  if (newstate.MX != MX_) {
1588  MVT::SetBlock(*newstate.MX,dimind,*MX_);
1589  }
1590  }
1591  }
1592  }
1593 
1594  // Compute R
1595  if(newstate.R == Teuchos::null)
1596  {
1597  om_->stream(Debug) << "Computing R\n";
1598 
1599  updateResidual();
1600  }
1601  // Copy R
1602  else
1603  {
1604  om_->stream(Debug) << "Copying R\n";
1605 
1606  if(computeAllRes_ == false)
1607  {
1608  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetGlobalLength(*newstate.R) != MVT::GetGlobalLength(*X_),
1609  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): vector length of newstate.R not correct." );
1610  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.R) != blockSize_,
1611  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): newstate.R must have at least block size vectors." );
1612  if (newstate.R != R_) {
1613  MVT::SetBlock(*newstate.R,bsind,*R_);
1614  }
1615  }
1616  else
1617  {
1618  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetGlobalLength(*newstate.R) != MVT::GetGlobalLength(*X_),
1619  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): vector length of newstate.R not correct." );
1620  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.R) != curDim_,
1621  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): newstate.R must have at least curDim vectors." );
1622  if (newstate.R != R_) {
1623  MVT::SetBlock(*newstate.R,dimind,*R_);
1624  }
1625  }
1626  }
1627 
1628  // R has been updated; mark the norms as out-of-date
1629  Rnorms_current_ = false;
1630  R2norms_current_ = false;
1631 
1632  // Set the largest safe shift
1633  largestSafeShift_ = newstate.largestSafeShift;
1634 
1635  // Copy over the Ritz shifts
1636  if(newstate.ritzShifts != Teuchos::null)
1637  {
1638  om_->stream(Debug) << "Copying Ritz shifts\n";
1639  std::copy(newstate.ritzShifts->begin(),newstate.ritzShifts->end(),ritzShifts_.begin());
1640  }
1641  else
1642  {
1643  om_->stream(Debug) << "Setting Ritz shifts to 0\n";
1644  for(size_t i=0; i<ritzShifts_.size(); i++)
1645  ritzShifts_[i] = ZERO;
1646  }
1647 
1648  for(size_t i=0; i<ritzShifts_.size(); i++)
1649  om_->stream(Debug) << "Ritz shifts[" << i << "] = " << ritzShifts_[i] << std::endl;
1650 
1651  // finally, we are initialized
1652  initialized_ = true;
1653 
1654  if (om_->isVerbosity( Debug ) ) {
1655  // Check almost everything here
1656  CheckList chk;
1657  chk.checkV = true;
1658  chk.checkX = true;
1659  chk.checkKX = true;
1660  chk.checkMX = true;
1661  chk.checkQ = true;
1662  chk.checkKK = true;
1663  om_->print( Debug, accuracyCheck(chk, ": after initialize()") );
1664  }
1665 
1666  // Print information on current status
1667  if (om_->isVerbosity(Debug)) {
1668  currentStatus( om_->stream(Debug) );
1669  }
1670  else if (om_->isVerbosity(IterationDetails)) {
1671  currentStatus( om_->stream(IterationDetails) );
1672  }
1673  }
1674 
1675 
1676 
1678  /* Initialize the state of the solver
1679  *
1680  * POST-CONDITIONS:
1681  *
1682  * V_ is orthonormal, orthogonal to auxVecs_, for first curDim_ vectors
1683  * theta_ contains Ritz w.r.t. V_(1:curDim_)
1684  * X is Ritz vectors w.r.t. V_(1:curDim_)
1685  * KX = Op*X
1686  * MX = M*X if hasM_
1687  * R = KX - MX*diag(theta_)
1688  *
1689  */
1690  template <class ScalarType, class MV, class OP>
1692  {
1693  // TODO: Check for bad input
1694  // NOTE: memory has been allocated by setBlockSize(). Use setBlock below; do not Clone
1695  // NOTE: Overall time spent in this routine is counted to timerInit_; portions will also be counted towards other primitives
1696 
1697  std::vector<int> bsind(blockSize_);
1698  for (int i=0; i<blockSize_; ++i) bsind[i] = i;
1699 
1700  // in TraceMinBase, V is primary
1701  // the order of dependence follows like so.
1702  // --init-> V,KK
1703  // --ritz analysis-> theta,X
1704  // --op apply-> KX,MX
1705  // --compute-> R
1706  //
1707  // if the user specifies all data for a level, we will accept it.
1708  // otherwise, we will generate the whole level, and all subsequent levels.
1709  //
1710  // the data members are ordered based on dependence, and the levels are
1711  // partitioned according to the amount of work required to produce the
1712  // items in a level.
1713  //
1714  // inconsistent multivectors widths and lengths will not be tolerated, and
1715  // will be treated with exceptions.
1716  //
1717  // for multivector pointers in newstate which point directly (as opposed to indirectly, via a view) to
1718  // multivectors in the solver, the copy will not be affected.
1719 
1720  // set up V and KK: get them from newstate if user specified them
1721  // otherwise, set them manually
1722  RCP<MV> lclV, lclKV, lclMV;
1724 
1725  // If V was supplied by the user...
1726  if (newstate.V != Teuchos::null) {
1727  om_->stream(Debug) << "Copying V from the user\n";
1728 
1729  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.V) != MVT::GetGlobalLength(*V_), std::invalid_argument,
1730  "Anasazi::TraceMinBase::initialize(newstate): Vector length of V not correct." );
1731  TEUCHOS_TEST_FOR_EXCEPTION( newstate.curDim < blockSize_, std::invalid_argument,
1732  "Anasazi::TraceMinBase::initialize(newstate): Rank of new state must be at least blockSize().");
1733  TEUCHOS_TEST_FOR_EXCEPTION( newstate.curDim > blockSize_*numBlocks_, std::invalid_argument,
1734  "Anasazi::TraceMinBase::initialize(newstate): Rank of new state must be less than getMaxSubspaceDim().");
1735 
1736  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.V) < newstate.curDim, std::invalid_argument,
1737  "Anasazi::TraceMinBase::initialize(newstate): Multivector for basis in new state must be as large as specified state rank.");
1738 
1739  curDim_ = newstate.curDim;
1740  // pick an integral amount
1741  curDim_ = (int)(curDim_ / blockSize_)*blockSize_;
1742 
1743  TEUCHOS_TEST_FOR_EXCEPTION( curDim_ != newstate.curDim, std::invalid_argument,
1744  "Anasazi::TraceMinBase::initialize(newstate): Rank of new state must be a multiple of getBlockSize().");
1745 
1746  // put data in V
1747  std::vector<int> nevind(curDim_);
1748  for (int i=0; i<curDim_; ++i) nevind[i] = i;
1749  if (newstate.V != V_) {
1750  if (curDim_ < MVT::GetNumberVecs(*newstate.V)) {
1751  newstate.V = MVT::CloneView(*newstate.V,nevind);
1752  }
1753  MVT::SetBlock(*newstate.V,nevind,*V_);
1754  }
1755  lclV = MVT::CloneViewNonConst(*V_,nevind);
1756  } // end if user specified V
1757  else
1758  {
1759  // get vectors from problem or generate something, projectAndNormalize
1760  RCP<const MV> ivec = problem_->getInitVec();
1761  TEUCHOS_TEST_FOR_EXCEPTION(ivec == Teuchos::null,std::invalid_argument,
1762  "Anasazi::TraceMinBase::initialize(newstate): Eigenproblem did not specify initial vectors to clone from.");
1763  // clear newstate so we won't use any data from it below
1764  newstate.X = Teuchos::null;
1765  newstate.MX = Teuchos::null;
1766  newstate.KX = Teuchos::null;
1767  newstate.R = Teuchos::null;
1768  newstate.T = Teuchos::null;
1769  newstate.RV = Teuchos::null;
1770  newstate.KK = Teuchos::null;
1771  newstate.KV = Teuchos::null;
1772  newstate.MopV = Teuchos::null;
1773  newstate.isOrtho = false;
1774 
1775  // If the user did not specify a current dimension, use that of the initial vectors they provided
1776  if(newstate.curDim > 0)
1777  curDim_ = newstate.curDim;
1778  else
1779  curDim_ = MVT::GetNumberVecs(*ivec);
1780 
1781  // pick the largest multiple of blockSize_
1782  curDim_ = (int)(curDim_ / blockSize_)*blockSize_;
1783  if (curDim_ > blockSize_*numBlocks_) {
1784  // user specified too many vectors... truncate
1785  // this produces a full subspace, but that is okay
1786  curDim_ = blockSize_*numBlocks_;
1787  }
1788  bool userand = false;
1789  if (curDim_ == 0) {
1790  // we need at least blockSize_ vectors
1791  // use a random multivec: ignore everything from InitVec
1792  userand = true;
1793  curDim_ = blockSize_;
1794  }
1795 
1796  std::vector<int> nevind(curDim_);
1797  for (int i=0; i<curDim_; ++i) nevind[i] = i;
1798 
1799  // get a pointer into V
1800  // lclV has curDim vectors
1801  //
1802  // get pointer to first curDim vectors in V_
1803  lclV = MVT::CloneViewNonConst(*V_,nevind);
1804  if (userand)
1805  {
1806  // generate random vector data
1807  MVT::MvRandom(*lclV);
1808  }
1809  else
1810  {
1811  if(newstate.curDim > 0)
1812  {
1813  if(MVT::GetNumberVecs(*newstate.V) > curDim_) {
1814  RCP<const MV> helperMV = MVT::CloneView(*newstate.V,nevind);
1815  MVT::SetBlock(*helperMV,nevind,*lclV);
1816  }
1817  // assign ivec to first part of lclV
1818  MVT::SetBlock(*newstate.V,nevind,*lclV);
1819  }
1820  else
1821  {
1822  if(MVT::GetNumberVecs(*ivec) > curDim_) {
1823  ivec = MVT::CloneView(*ivec,nevind);
1824  }
1825  // assign ivec to first part of lclV
1826  MVT::SetBlock(*ivec,nevind,*lclV);
1827  }
1828  }
1829  } // end if user did not specify V
1830 
1831  // Nuke everything from orbit
1832  // This is a temporary measure due to a bug in the code that I have not found yet
1833  // It adds a minimal amount of work
1834  newstate.X = Teuchos::null;
1835  newstate.MX = Teuchos::null;
1836  newstate.KX = Teuchos::null;
1837  newstate.R = Teuchos::null;
1838  newstate.T = Teuchos::null;
1839  newstate.RV = Teuchos::null;
1840  newstate.KK = Teuchos::null;
1841  newstate.KV = Teuchos::null;
1842  newstate.MopV = Teuchos::null;
1843  newstate.isOrtho = false;
1844 
1845  // A vector of relevant indices
1846  std::vector<int> dimind(curDim_);
1847  for (int i=0; i<curDim_; ++i) dimind[i] = i;
1848 
1849  // Project the auxVecs out of V
1850  if(auxVecs_.size() > 0)
1851  orthman_->projectMat(*lclV,auxVecs_);
1852 
1853  // Compute KV
1854  if(Op_ != Teuchos::null && newstate.KV == Teuchos::null)
1855  {
1856  om_->stream(Debug) << "Computing KV\n";
1857 
1858 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1859  Teuchos::TimeMonitor lcltimer( *timerOp_ );
1860 #endif
1861  count_ApplyOp_+= curDim_;
1862 
1863  // These things are now invalid
1864  newstate.X = Teuchos::null;
1865  newstate.MX = Teuchos::null;
1866  newstate.KX = Teuchos::null;
1867  newstate.R = Teuchos::null;
1868  newstate.T = Teuchos::null;
1869  newstate.RV = Teuchos::null;
1870  newstate.KK = Teuchos::null;
1871 
1872  lclKV = MVT::CloneViewNonConst(*KV_,dimind);
1873  OPT::Apply(*Op_,*lclV,*lclKV);
1874  }
1875  // Copy KV
1876  else if(Op_ != Teuchos::null)
1877  {
1878  om_->stream(Debug) << "Copying KV\n";
1879 
1880  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.KV) != MVT::GetGlobalLength(*KV_), std::invalid_argument,
1881  "Anasazi::TraceMinBase::initialize(newstate): Vector length of KV not correct." );
1882 
1883  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.KV) < curDim_, std::invalid_argument,
1884  "Anasazi::TraceMinBase::initialize(newstate): Number of vectors in KV not correct.");
1885 
1886  if (newstate.KV != KV_) {
1887  if (curDim_ < MVT::GetNumberVecs(*newstate.KV)) {
1888  newstate.KV = MVT::CloneView(*newstate.KV,dimind);
1889  }
1890  MVT::SetBlock(*newstate.KV,dimind,*KV_);
1891  }
1892  lclKV = MVT::CloneViewNonConst(*KV_,dimind);
1893  }
1894  else
1895  {
1896  om_->stream(Debug) << "There is no KV\n";
1897 
1898  lclKV = lclV;
1899  KV_ = V_;
1900  }
1901 
1902 
1903 
1904  // Project and normalize so that V'V = I and Q'V=0
1905  if(newstate.isOrtho == false)
1906  {
1907  om_->stream(Debug) << "Project and normalize\n";
1908 
1909 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1910  Teuchos::TimeMonitor lcltimer( *timerOrtho_ );
1911 #endif
1912 
1913  // These things are now invalid
1914  newstate.MopV = Teuchos::null;
1915  newstate.X = Teuchos::null;
1916  newstate.MX = Teuchos::null;
1917  newstate.KX = Teuchos::null;
1918  newstate.R = Teuchos::null;
1919  newstate.T = Teuchos::null;
1920  newstate.RV = Teuchos::null;
1921  newstate.KK = Teuchos::null;
1922 
1923  // Normalize lclKV
1925  int rank = orthman_->normalizeMat(*lclKV,gamma);
1926 
1927  // lclV = lclV/gamma
1929  SDsolver.setMatrix(gamma);
1930  SDsolver.invert();
1931  RCP<MV> tempMV = MVT::CloneCopy(*lclV);
1932  MVT::MvTimesMatAddMv(ONE,*tempMV,*gamma,ZERO,*lclV);
1933 
1934  TEUCHOS_TEST_FOR_EXCEPTION(rank != curDim_,TraceMinBaseInitFailure,
1935  "Anasazi::TraceMinBase::initialize(): Couldn't generate initial basis of full rank.");
1936  }
1937 
1938  // Compute MV if necessary
1939  if(hasM_ && newstate.MopV == Teuchos::null)
1940  {
1941  om_->stream(Debug) << "Computing MV\n";
1942 
1943 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1944  Teuchos::TimeMonitor lcltimer( *timerMOp_ );
1945 #endif
1946  count_ApplyM_+= curDim_;
1947 
1948  // Get a pointer to the relevant parts of MV
1949  lclMV = MVT::CloneViewNonConst(*MV_,dimind);
1950  OPT::Apply(*MOp_,*lclV,*lclMV);
1951  }
1952  // Copy MV if necessary
1953  else if(hasM_)
1954  {
1955  om_->stream(Debug) << "Copying MV\n";
1956 
1957  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.MopV) != MVT::GetGlobalLength(*MV_), std::invalid_argument,
1958  "Anasazi::TraceMinBase::initialize(newstate): Vector length of MV not correct." );
1959 
1960  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.MopV) < curDim_, std::invalid_argument,
1961  "Anasazi::TraceMinBase::initialize(newstate): Number of vectors in MV not correct.");
1962 
1963  if(newstate.MopV != MV_) {
1964  if(curDim_ < MVT::GetNumberVecs(*newstate.MopV)) {
1965  newstate.MopV = MVT::CloneView(*newstate.MopV,dimind);
1966  }
1967  MVT::SetBlock(*newstate.MopV,dimind,*MV_);
1968  }
1969  lclMV = MVT::CloneViewNonConst(*MV_,dimind);
1970  }
1971  // There is no M, so there is no MV
1972  else
1973  {
1974  om_->stream(Debug) << "There is no MV\n";
1975 
1976  lclMV = lclV;
1977  MV_ = V_;
1978  }
1979 
1980  // Compute KK
1981  if(newstate.KK == Teuchos::null)
1982  {
1983  om_->stream(Debug) << "Computing KK\n";
1984 
1985  // These things are now invalid
1986  newstate.X = Teuchos::null;
1987  newstate.MX = Teuchos::null;
1988  newstate.KX = Teuchos::null;
1989  newstate.R = Teuchos::null;
1990  newstate.T = Teuchos::null;
1991  newstate.RV = Teuchos::null;
1992 
1993  // Note: there used to be a bug here.
1994  // We can't just call computeKK because it only computes the new parts of KK
1995 
1996  // Get a pointer to the part of KK we're interested in
1997  lclKK = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,curDim_) );
1998 
1999  // KK := V'KV
2000  MVT::MvTransMv(ONE,*lclV,*lclKV,*lclKK);
2001  }
2002  // Copy KK
2003  else
2004  {
2005  om_->stream(Debug) << "Copying KK\n";
2006 
2007  // check size of KK
2008  TEUCHOS_TEST_FOR_EXCEPTION( newstate.KK->numRows() < curDim_ || newstate.KK->numCols() < curDim_, std::invalid_argument,
2009  "Anasazi::TraceMinBase::initialize(newstate): Projected matrix in new state must be as large as specified state rank.");
2010 
2011  // put data into KK_
2012  lclKK = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,curDim_) );
2013  if (newstate.KK != KK_) {
2014  if (newstate.KK->numRows() > curDim_ || newstate.KK->numCols() > curDim_) {
2015  newstate.KK = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*newstate.KK,curDim_,curDim_) );
2016  }
2017  lclKK->assign(*newstate.KK);
2018  }
2019  }
2020 
2021  // Compute Ritz pairs
2022  if(newstate.T == Teuchos::null || newstate.RV == Teuchos::null)
2023  {
2024  om_->stream(Debug) << "Computing Ritz pairs\n";
2025 
2026  // These things are now invalid
2027  newstate.X = Teuchos::null;
2028  newstate.MX = Teuchos::null;
2029  newstate.KX = Teuchos::null;
2030  newstate.R = Teuchos::null;
2031  newstate.T = Teuchos::null;
2032  newstate.RV = Teuchos::null;
2033 
2034  computeRitzPairs();
2035  }
2036  // Copy Ritz pairs
2037  else
2038  {
2039  om_->stream(Debug) << "Copying Ritz pairs\n";
2040 
2041  TEUCHOS_TEST_FOR_EXCEPTION((signed int)(newstate.T->size()) != curDim_,
2042  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of T must be consistent with dimension of V.");
2043 
2044  TEUCHOS_TEST_FOR_EXCEPTION( newstate.RV->numRows() < curDim_ || newstate.RV->numCols() < curDim_, std::invalid_argument,
2045  "Anasazi::TraceMinBase::initialize(newstate): Ritz vectors in new state must be as large as specified state rank.");
2046 
2047  std::copy(newstate.T->begin(),newstate.T->end(),theta_.begin());
2048 
2049  lclRV = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,curDim_) );
2050  if (newstate.RV != ritzVecs_) {
2051  if (newstate.RV->numRows() > curDim_ || newstate.RV->numCols() > curDim_) {
2052  newstate.RV = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*newstate.RV,curDim_,curDim_) );
2053  }
2054  lclRV->assign(*newstate.RV);
2055  }
2056  }
2057 
2058  // Compute X
2059  if(newstate.X == Teuchos::null)
2060  {
2061  om_->stream(Debug) << "Computing X\n";
2062 
2063  // These things are now invalid
2064  newstate.MX = Teuchos::null;
2065  newstate.KX = Teuchos::null;
2066  newstate.R = Teuchos::null;
2067 
2068  computeX();
2069  }
2070  // Copy X
2071  else
2072  {
2073  om_->stream(Debug) << "Copying X\n";
2074 
2075  if(computeAllRes_ == false)
2076  {
2077  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.X) != blockSize_ || MVT::GetGlobalLength(*newstate.X) != MVT::GetGlobalLength(*X_),
2078  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of X must be consistent with block size and length of V.");
2079 
2080  if(newstate.X != X_) {
2081  MVT::SetBlock(*newstate.X,bsind,*X_);
2082  }
2083  }
2084  else
2085  {
2086  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.X) != curDim_ || MVT::GetGlobalLength(*newstate.X) != MVT::GetGlobalLength(*X_),
2087  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of X must be consistent with current dimension and length of V.");
2088 
2089  if(newstate.X != X_) {
2090  MVT::SetBlock(*newstate.X,dimind,*X_);
2091  }
2092  }
2093  }
2094 
2095  // Compute KX and MX if necessary
2096  // TODO: These technically should be separate; it won't matter much in terms of running time though
2097  if((Op_ != Teuchos::null && newstate.KX == Teuchos::null) || (hasM_ && newstate.MX == Teuchos::null))
2098  {
2099  om_->stream(Debug) << "Computing KX and MX\n";
2100 
2101  // These things are now invalid
2102  newstate.R = Teuchos::null;
2103 
2104  updateKXMX();
2105  }
2106  // Copy KX and MX if necessary
2107  else
2108  {
2109  om_->stream(Debug) << "Copying KX and MX\n";
2110  if(Op_ != Teuchos::null)
2111  {
2112  if(computeAllRes_ == false)
2113  {
2114  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.KX) != blockSize_ || MVT::GetGlobalLength(*newstate.KX) != MVT::GetGlobalLength(*KX_),
2115  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of KX must be consistent with block size and length of V.");
2116 
2117  if(newstate.KX != KX_) {
2118  MVT::SetBlock(*newstate.KX,bsind,*KX_);
2119  }
2120  }
2121  else
2122  {
2123  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.KX) != curDim_ || MVT::GetGlobalLength(*newstate.KX) != MVT::GetGlobalLength(*KX_),
2124  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of KX must be consistent with current dimension and length of V.");
2125 
2126  if (newstate.KX != KX_) {
2127  MVT::SetBlock(*newstate.KX,dimind,*KX_);
2128  }
2129  }
2130  }
2131 
2132  if(hasM_)
2133  {
2134  if(computeAllRes_ == false)
2135  {
2136  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.MX) != blockSize_ || MVT::GetGlobalLength(*newstate.MX) != MVT::GetGlobalLength(*MX_),
2137  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of MX must be consistent with block size and length of V.");
2138 
2139  if (newstate.MX != MX_) {
2140  MVT::SetBlock(*newstate.MX,bsind,*MX_);
2141  }
2142  }
2143  else
2144  {
2145  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.MX) != curDim_ || MVT::GetGlobalLength(*newstate.MX) != MVT::GetGlobalLength(*MX_),
2146  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): Size of MX must be consistent with current dimension and length of V.");
2147 
2148  if (newstate.MX != MX_) {
2149  MVT::SetBlock(*newstate.MX,dimind,*MX_);
2150  }
2151  }
2152  }
2153  }
2154 
2155  // Scale X so each vector is of length 1
2156  {
2157  // Get norm of each vector in X
2158  const int nvecs = computeAllRes_ ? curDim_ : blockSize_;
2159  Teuchos::Range1D dimind2 (0, nvecs-1);
2160  RCP<MV> lclX = MVT::CloneViewNonConst(*X_, dimind2);
2161  std::vector<ScalarType> normvec(nvecs);
2162  orthman_->normMat(*lclX,normvec);
2163 
2164  // Scale X, KX, and MX accordingly
2165  for (int i = 0; i < nvecs; ++i) {
2166  normvec[i] = ONE / normvec[i];
2167  }
2168  MVT::MvScale (*lclX, normvec);
2169  if (Op_ != Teuchos::null) {
2170  RCP<MV> lclKX = MVT::CloneViewNonConst (*KX_, dimind2);
2171  MVT::MvScale (*lclKX, normvec);
2172  }
2173  if (hasM_) {
2174  RCP<MV> lclMX = MVT::CloneViewNonConst (*MX_, dimind2);
2175  MVT::MvScale (*lclMX, normvec);
2176  }
2177 
2178  // Scale eigenvalues
2179  for (int i = 0; i < nvecs; ++i) {
2180  theta_[i] = theta_[i] * normvec[i] * normvec[i];
2181  }
2182  }
2183 
2184  // Compute R
2185  if(newstate.R == Teuchos::null)
2186  {
2187  om_->stream(Debug) << "Computing R\n";
2188 
2189  updateResidual();
2190  }
2191  // Copy R
2192  else
2193  {
2194  om_->stream(Debug) << "Copying R\n";
2195 
2196  if(computeAllRes_ == false)
2197  {
2198  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetGlobalLength(*newstate.R) != MVT::GetGlobalLength(*X_),
2199  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): vector length of newstate.R not correct." );
2200  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.R) != blockSize_,
2201  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): newstate.R must have at least block size vectors." );
2202  if (newstate.R != R_) {
2203  MVT::SetBlock(*newstate.R,bsind,*R_);
2204  }
2205  }
2206  else
2207  {
2208  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetGlobalLength(*newstate.R) != MVT::GetGlobalLength(*X_),
2209  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): vector length of newstate.R not correct." );
2210  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.R) != curDim_,
2211  std::invalid_argument, "Anasazi::TraceMinBase::initialize(newstate): newstate.R must have at least curDim vectors." );
2212  if (newstate.R != R_) {
2213  MVT::SetBlock(*newstate.R,dimind,*R_);
2214  }
2215  }
2216  }
2217 
2218  // R has been updated; mark the norms as out-of-date
2219  Rnorms_current_ = false;
2220  R2norms_current_ = false;
2221 
2222  // Set the largest safe shift
2223  largestSafeShift_ = newstate.largestSafeShift;
2224 
2225  // Copy over the Ritz shifts
2226  if(newstate.ritzShifts != Teuchos::null)
2227  {
2228  om_->stream(Debug) << "Copying Ritz shifts\n";
2229  std::copy(newstate.ritzShifts->begin(),newstate.ritzShifts->end(),ritzShifts_.begin());
2230  }
2231  else
2232  {
2233  om_->stream(Debug) << "Setting Ritz shifts to 0\n";
2234  for(size_t i=0; i<ritzShifts_.size(); i++)
2235  ritzShifts_[i] = ZERO;
2236  }
2237 
2238  for(size_t i=0; i<ritzShifts_.size(); i++)
2239  om_->stream(Debug) << "Ritz shifts[" << i << "] = " << ritzShifts_[i] << std::endl;
2240 
2241  // finally, we are initialized
2242  initialized_ = true;
2243 
2244  if (om_->isVerbosity( Debug ) ) {
2245  // Check almost everything here
2246  CheckList chk;
2247  chk.checkV = true;
2248  chk.checkX = true;
2249  chk.checkKX = true;
2250  chk.checkMX = true;
2251  chk.checkQ = true;
2252  chk.checkKK = true;
2253  om_->print( Debug, accuracyCheck(chk, ": after initialize()") );
2254  }
2255 
2256  // Print information on current status
2257  if (om_->isVerbosity(Debug)) {
2258  currentStatus( om_->stream(Debug) );
2259  }
2260  else if (om_->isVerbosity(IterationDetails)) {
2261  currentStatus( om_->stream(IterationDetails) );
2262  }
2263  }
2264 
2265 
2267  // Return initialized state
2268  template <class ScalarType, class MV, class OP>
2269  bool TraceMinBase<ScalarType,MV,OP>::isInitialized() const { return initialized_; }
2270 
2271 
2273  // Set the block size and make necessary adjustments.
2274  template <class ScalarType, class MV, class OP>
2275  void TraceMinBase<ScalarType,MV,OP>::setSize (int blockSize, int numBlocks)
2276  {
2277  // This routine only allocates space; it doesn't not perform any computation
2278  // any change in size will invalidate the state of the solver.
2279 
2280  TEUCHOS_TEST_FOR_EXCEPTION(blockSize < 1, std::invalid_argument, "Anasazi::TraceMinBase::setSize(blocksize,numblocks): blocksize must be strictly positive.");
2281 
2282  if (blockSize == blockSize_ && numBlocks == numBlocks_) {
2283  // do nothing
2284  return;
2285  }
2286 
2287  blockSize_ = blockSize;
2288  numBlocks_ = numBlocks;
2289 
2290  RCP<const MV> tmp;
2291  // grab some Multivector to Clone
2292  // in practice, getInitVec() should always provide this, but it is possible to use a
2293  // Eigenproblem with nothing in getInitVec() by manually initializing with initialize();
2294  // in case of that strange scenario, we will try to Clone from X_ first, then resort to getInitVec()
2295  if (X_ != Teuchos::null) { // this is equivalent to blockSize_ > 0
2296  tmp = X_;
2297  }
2298  else {
2299  tmp = problem_->getInitVec();
2300  TEUCHOS_TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::invalid_argument,
2301  "Anasazi::TraceMinBase::setSize(): eigenproblem did not specify initial vectors to clone from.");
2302  }
2303 
2304  TEUCHOS_TEST_FOR_EXCEPTION(numAuxVecs_+blockSize*static_cast<ptrdiff_t>(numBlocks) > MVT::GetGlobalLength(*tmp),std::invalid_argument,
2305  "Anasazi::TraceMinBase::setSize(): max subspace dimension and auxilliary subspace too large. Potentially impossible orthogonality constraints.");
2306 
2307  // New subspace dimension
2308  int newsd = blockSize_*numBlocks_;
2309 
2311  // blockSize dependent
2312  //
2313  ritzShifts_.resize(blockSize_,ZERO);
2314  if(computeAllRes_ == false)
2315  {
2316  // grow/allocate vectors
2317  Rnorms_.resize(blockSize_,NANVAL);
2318  R2norms_.resize(blockSize_,NANVAL);
2319  //
2320  // clone multivectors off of tmp
2321  //
2322  // free current allocation first, to make room for new allocation
2323  X_ = Teuchos::null;
2324  KX_ = Teuchos::null;
2325  MX_ = Teuchos::null;
2326  R_ = Teuchos::null;
2327  V_ = Teuchos::null;
2328  KV_ = Teuchos::null;
2329  MV_ = Teuchos::null;
2330 
2331  om_->print(Debug," >> Allocating X_\n");
2332  X_ = MVT::Clone(*tmp,blockSize_);
2333  if(Op_ != Teuchos::null) {
2334  om_->print(Debug," >> Allocating KX_\n");
2335  KX_ = MVT::Clone(*tmp,blockSize_);
2336  }
2337  else {
2338  KX_ = X_;
2339  }
2340  if (hasM_) {
2341  om_->print(Debug," >> Allocating MX_\n");
2342  MX_ = MVT::Clone(*tmp,blockSize_);
2343  }
2344  else {
2345  MX_ = X_;
2346  }
2347  om_->print(Debug," >> Allocating R_\n");
2348  R_ = MVT::Clone(*tmp,blockSize_);
2349  }
2350  else
2351  {
2352  // grow/allocate vectors
2353  Rnorms_.resize(newsd,NANVAL);
2354  R2norms_.resize(newsd,NANVAL);
2355  //
2356  // clone multivectors off of tmp
2357  //
2358  // free current allocation first, to make room for new allocation
2359  X_ = Teuchos::null;
2360  KX_ = Teuchos::null;
2361  MX_ = Teuchos::null;
2362  R_ = Teuchos::null;
2363  V_ = Teuchos::null;
2364  KV_ = Teuchos::null;
2365  MV_ = Teuchos::null;
2366 
2367  om_->print(Debug," >> Allocating X_\n");
2368  X_ = MVT::Clone(*tmp,newsd);
2369  if(Op_ != Teuchos::null) {
2370  om_->print(Debug," >> Allocating KX_\n");
2371  KX_ = MVT::Clone(*tmp,newsd);
2372  }
2373  else {
2374  KX_ = X_;
2375  }
2376  if (hasM_) {
2377  om_->print(Debug," >> Allocating MX_\n");
2378  MX_ = MVT::Clone(*tmp,newsd);
2379  }
2380  else {
2381  MX_ = X_;
2382  }
2383  om_->print(Debug," >> Allocating R_\n");
2384  R_ = MVT::Clone(*tmp,newsd);
2385  }
2386 
2388  // blockSize*numBlocks dependent
2389  //
2390  theta_.resize(newsd,NANVAL);
2391  om_->print(Debug," >> Allocating V_\n");
2392  V_ = MVT::Clone(*tmp,newsd);
2393  KK_ = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(newsd,newsd) );
2394  ritzVecs_ = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(newsd,newsd) );
2395 
2396  if(Op_ != Teuchos::null) {
2397  om_->print(Debug," >> Allocating KV_\n");
2398  KV_ = MVT::Clone(*tmp,newsd);
2399  }
2400  else {
2401  KV_ = V_;
2402  }
2403  if (hasM_) {
2404  om_->print(Debug," >> Allocating MV_\n");
2405  MV_ = MVT::Clone(*tmp,newsd);
2406  }
2407  else {
2408  MV_ = V_;
2409  }
2410 
2411  om_->print(Debug," >> done allocating.\n");
2412 
2413  initialized_ = false;
2414  curDim_ = 0;
2415  }
2416 
2417 
2419  // Set the auxiliary vectors
2420  template <class ScalarType, class MV, class OP>
2422  typedef typename Teuchos::Array<RCP<const MV> >::iterator tarcpmv;
2423 
2424  // set new auxiliary vectors
2425  auxVecs_ = auxvecs;
2426 
2427  if(hasM_)
2428  MauxVecs_.resize(0);
2429  numAuxVecs_ = 0;
2430 
2431  for (tarcpmv i=auxVecs_.begin(); i != auxVecs_.end(); ++i) {
2432  numAuxVecs_ += MVT::GetNumberVecs(**i);
2433 
2434  if(hasM_)
2435  {
2436 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
2437  Teuchos::TimeMonitor lcltimer( *timerMOp_ );
2438 #endif
2439  count_ApplyM_+= MVT::GetNumberVecs(**i);
2440 
2441  RCP<MV> helperMV = MVT::Clone(**i,MVT::GetNumberVecs(**i));
2442  OPT::Apply(*MOp_,**i,*helperMV);
2443  MauxVecs_.push_back(helperMV);
2444  }
2445  }
2446 
2447  // If the solver has been initialized, V is not necessarily orthogonal to new auxiliary vectors
2448  if (numAuxVecs_ > 0 && initialized_) {
2449  initialized_ = false;
2450  }
2451 
2452  if (om_->isVerbosity( Debug ) ) {
2453  // Check almost everything here
2454  CheckList chk;
2455  chk.checkQ = true;
2456  om_->print( Debug, accuracyCheck(chk, ": after setAuxVecs()") );
2457  }
2458  }
2459 
2460 
2462  // compute/return residual M-norms
2463  template <class ScalarType, class MV, class OP>
2464  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
2466  if (Rnorms_current_ == false) {
2467  // Update the residual norms
2468  if(computeAllRes_)
2469  {
2470  std::vector<int> curind(curDim_);
2471  for(int i=0; i<curDim_; i++)
2472  curind[i] = i;
2473 
2474  RCP<const MV> locR = MVT::CloneView(*R_,curind);
2475  std::vector<ScalarType> locNorms(curDim_);
2476  orthman_->norm(*locR,locNorms);
2477 
2478  for(int i=0; i<curDim_; i++)
2479  Rnorms_[i] = locNorms[i];
2480  for(int i=curDim_+1; i<blockSize_*numBlocks_; i++)
2481  Rnorms_[i] = NANVAL;
2482 
2483  Rnorms_current_ = true;
2484  locNorms.resize(blockSize_);
2485  return locNorms;
2486  }
2487  else
2488  orthman_->norm(*R_,Rnorms_);
2489  Rnorms_current_ = true;
2490  }
2491  else if(computeAllRes_)
2492  {
2493  std::vector<ScalarType> locNorms(blockSize_);
2494  for(int i=0; i<blockSize_; i++)
2495  locNorms[i] = Rnorms_[i];
2496  return locNorms;
2497  }
2498 
2499  return Rnorms_;
2500  }
2501 
2502 
2504  // compute/return residual 2-norms
2505  template <class ScalarType, class MV, class OP>
2506  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
2508  if (R2norms_current_ == false) {
2509  // Update the residual 2-norms
2510  if(computeAllRes_)
2511  {
2512  std::vector<int> curind(curDim_);
2513  for(int i=0; i<curDim_; i++)
2514  curind[i] = i;
2515 
2516  RCP<const MV> locR = MVT::CloneView(*R_,curind);
2517  std::vector<ScalarType> locNorms(curDim_);
2518  MVT::MvNorm(*locR,locNorms);
2519 
2520  for(int i=0; i<curDim_; i++)
2521  {
2522  R2norms_[i] = locNorms[i];
2523  }
2524  for(int i=curDim_+1; i<blockSize_*numBlocks_; i++)
2525  R2norms_[i] = NANVAL;
2526 
2527  R2norms_current_ = true;
2528  locNorms.resize(blockSize_);
2529  return locNorms;
2530  }
2531  else
2532  MVT::MvNorm(*R_,R2norms_);
2533  R2norms_current_ = true;
2534  }
2535  else if(computeAllRes_)
2536  {
2537  std::vector<ScalarType> locNorms(blockSize_);
2538  for(int i=0; i<blockSize_; i++)
2539  locNorms[i] = R2norms_[i];
2540  return locNorms;
2541  }
2542 
2543  return R2norms_;
2544  }
2545 
2547  // Set a new StatusTest for the solver.
2548  template <class ScalarType, class MV, class OP>
2550  TEUCHOS_TEST_FOR_EXCEPTION(test == Teuchos::null,std::invalid_argument,
2551  "Anasazi::TraceMinBase::setStatusTest() was passed a null StatusTest.");
2552  tester_ = test;
2553  }
2554 
2556  // Get the current StatusTest used by the solver.
2557  template <class ScalarType, class MV, class OP>
2559  return tester_;
2560  }
2561 
2563  // Print the current status of the solver
2564  template <class ScalarType, class MV, class OP>
2565  void
2567  {
2568  using std::endl;
2569 
2570  os.setf(std::ios::scientific, std::ios::floatfield);
2571  os.precision(6);
2572  os <<endl;
2573  os <<"================================================================================" << endl;
2574  os << endl;
2575  os <<" TraceMinBase Solver Status" << endl;
2576  os << endl;
2577  os <<"The solver is "<<(initialized_ ? "initialized." : "not initialized.") << endl;
2578  os <<"The number of iterations performed is " <<iter_<<endl;
2579  os <<"The block size is " << blockSize_<<endl;
2580  os <<"The number of blocks is " << numBlocks_<<endl;
2581  os <<"The current basis size is " << curDim_<<endl;
2582  os <<"The number of auxiliary vectors is "<< numAuxVecs_ << endl;
2583  os <<"The number of operations Op*x is "<<count_ApplyOp_<<endl;
2584  os <<"The number of operations M*x is "<<count_ApplyM_<<endl;
2585 
2586  os.setf(std::ios_base::right, std::ios_base::adjustfield);
2587 
2588  if (initialized_) {
2589  os << endl;
2590  os <<"CURRENT EIGENVALUE ESTIMATES "<<endl;
2591  os << std::setw(20) << "Eigenvalue"
2592  << std::setw(20) << "Residual(M)"
2593  << std::setw(20) << "Residual(2)"
2594  << endl;
2595  os <<"--------------------------------------------------------------------------------"<<endl;
2596  for (int i=0; i<blockSize_; ++i) {
2597  os << std::setw(20) << theta_[i];
2598  if (Rnorms_current_) os << std::setw(20) << Rnorms_[i];
2599  else os << std::setw(20) << "not current";
2600  if (R2norms_current_) os << std::setw(20) << R2norms_[i];
2601  else os << std::setw(20) << "not current";
2602  os << endl;
2603  }
2604  }
2605  os <<"================================================================================" << endl;
2606  os << endl;
2607  }
2608 
2610  template <class ScalarType, class MV, class OP>
2611  ScalarType TraceMinBase<ScalarType,MV,OP>::getTrace() const
2612  {
2613  ScalarType currentTrace = ZERO;
2614 
2615  for(int i=0; i < blockSize_; i++)
2616  currentTrace += theta_[i];
2617 
2618  return currentTrace;
2619  }
2620 
2622  template <class ScalarType, class MV, class OP>
2623  bool TraceMinBase<ScalarType,MV,OP>::traceLeveled()
2624  {
2625  ScalarType ratioOfChange = traceThresh_;
2626 
2627  if(previouslyLeveled_)
2628  {
2629  om_->stream(Debug) << "The trace already leveled, so we're not going to check it again\n";
2630  return true;
2631  }
2632 
2633  ScalarType currentTrace = getTrace();
2634 
2635  om_->stream(Debug) << "The current trace is " << currentTrace << std::endl;
2636 
2637  // Compute the ratio of the change
2638  // We seek the point where the trace has leveled off
2639  // It should be reasonably safe to shift at this point
2640  if(previousTrace_ != ZERO)
2641  {
2642  om_->stream(Debug) << "The previous trace was " << previousTrace_ << std::endl;
2643 
2644  ratioOfChange = std::abs(previousTrace_-currentTrace)/std::abs(previousTrace_);
2645  om_->stream(Debug) << "The ratio of change is " << ratioOfChange << std::endl;
2646  }
2647 
2648  previousTrace_ = currentTrace;
2649 
2650  if(ratioOfChange < traceThresh_)
2651  {
2652  previouslyLeveled_ = true;
2653  return true;
2654  }
2655  return false;
2656  }
2657 
2659  // Compute the residual of each CLUSTER of eigenvalues
2660  // This is important for selecting the Ritz shifts
2661  template <class ScalarType, class MV, class OP>
2662  std::vector<ScalarType> TraceMinBase<ScalarType,MV,OP>::getClusterResids()
2663  {
2664  int nvecs;
2665 
2666  if(computeAllRes_)
2667  nvecs = curDim_;
2668  else
2669  nvecs = blockSize_;
2670 
2671  getRes2Norms();
2672 
2673  std::vector<ScalarType> clusterResids(nvecs);
2674  std::vector<int> clusterIndices;
2675  if(considerClusters_)
2676  {
2677  for(int i=0; i < nvecs; i++)
2678  {
2679  // test for cluster
2680  if(clusterIndices.empty() || (theta_[i-1] + R2norms_[i-1] >= theta_[i] - R2norms_[i]))
2681  {
2682  // Add to cluster
2683  if(!clusterIndices.empty()) om_->stream(Debug) << theta_[i-1] << " is in a cluster with " << theta_[i] << " because " << theta_[i-1] + R2norms_[i-1] << " >= " << theta_[i] - R2norms_[i] << std::endl;
2684  clusterIndices.push_back(i);
2685  }
2686  // Cluster completed
2687  else
2688  {
2689  om_->stream(Debug) << theta_[i-1] << " is NOT in a cluster with " << theta_[i] << " because " << theta_[i-1] + R2norms_[i-1] << " < " << theta_[i] - R2norms_[i] << std::endl;
2690  ScalarType totalRes = ZERO;
2691  for(size_t j=0; j < clusterIndices.size(); j++)
2692  totalRes += (R2norms_[clusterIndices[j]]*R2norms_[clusterIndices[j]]);
2693 
2694  // If the smallest magnitude value of this sign is in a cluster with the
2695  // largest magnitude cluster of this sign, it is not safe for the smallest
2696  // eigenvalue to use a shift
2697  if(theta_[clusterIndices[0]] < 0 && theta_[i] < 0)
2698  negSafeToShift_ = true;
2699  else if(theta_[clusterIndices[0]] > 0 && theta_[i] > 0)
2700  posSafeToShift_ = true;
2701 
2702  for(size_t j=0; j < clusterIndices.size(); j++)
2703  clusterResids[clusterIndices[j]] = sqrt(totalRes);
2704 
2705  clusterIndices.clear();
2706  clusterIndices.push_back(i);
2707  }
2708  }
2709 
2710  // Handle last cluster
2711  ScalarType totalRes = ZERO;
2712  for(size_t j=0; j < clusterIndices.size(); j++)
2713  totalRes += R2norms_[clusterIndices[j]];
2714  for(size_t j=0; j < clusterIndices.size(); j++)
2715  clusterResids[clusterIndices[j]] = totalRes;
2716  }
2717  else
2718  {
2719  for(int j=0; j < nvecs; j++)
2720  clusterResids[j] = R2norms_[j];
2721  }
2722 
2723  return clusterResids;
2724  }
2725 
2727  // Compute the Ritz shifts based on the Ritz values and residuals
2728  // A note on shifting: if the matrix is indefinite, you NEED to use a large block size
2729  // TODO: resids[i] on its own is unsafe for the generalized EVP
2730  // See "A Parallel Implementation of the Trace Minimization Eigensolver"
2731  // by Eloy Romero and Jose E. Roman
2732  template <class ScalarType, class MV, class OP>
2733  void TraceMinBase<ScalarType,MV,OP>::computeRitzShifts(const std::vector<ScalarType>& clusterResids)
2734  {
2735  std::vector<ScalarType> thetaMag(theta_);
2736  bool traceHasLeveled = traceLeveled();
2737 
2738  // Compute the magnitude of the eigenvalues
2739  for(int i=0; i<blockSize_; i++)
2740  thetaMag[i] = std::abs(thetaMag[i]);
2741 
2742  // Test whether it is safe to shift
2743  // TODO: Add residual shift option
2744  // Note: If you choose single shift with an indefinite matrix, you're gonna have a bad time...
2745  // Note: This is not safe for indefinite matrices, and I don't even know that it CAN be made safe.
2746  // There just isn't any theoretical reason it should work.
2747  // TODO: It feels like this should be different for TraceMin-Base; we should be able to use all eigenvalues in the current subspace to determine whether we have a cluster
2748  if(whenToShift_ == ALWAYS_SHIFT || (whenToShift_ == SHIFT_WHEN_TRACE_LEVELS && traceHasLeveled))
2749  {
2750  // Set the shift to the largest safe shift
2751  if(howToShift_ == LARGEST_CONVERGED_SHIFT)
2752  {
2753  for(int i=0; i<blockSize_; i++)
2754  ritzShifts_[i] = largestSafeShift_;
2755  }
2756  // Set the shifts to the Ritz values
2757  else if(howToShift_ == RITZ_VALUES_SHIFT)
2758  {
2759  ritzShifts_[0] = theta_[0];
2760 
2761  // If we're using mulitple shifts, set them to EACH Ritz value.
2762  // Otherwise, only use the smallest
2763  if(useMultipleShifts_)
2764  {
2765  for(int i=1; i<blockSize_; i++)
2766  ritzShifts_[i] = theta_[i];
2767  }
2768  else
2769  {
2770  for(int i=1; i<blockSize_; i++)
2771  ritzShifts_[i] = ritzShifts_[0];
2772  }
2773  }
2774  else if(howToShift_ == EXPERIMENTAL_SHIFT)
2775  {
2776  ritzShifts_[0] = std::max(largestSafeShift_,theta_[0]-clusterResids[0]);
2777  for(int i=1; i<blockSize_; i++)
2778  {
2779  ritzShifts_[i] = std::max(ritzShifts_[i-1],theta_[i]-clusterResids[i]);
2780  }
2781  }
2782  // Use Dr. Sameh's original shifting strategy
2783  else if(howToShift_ == ADJUSTED_RITZ_SHIFT)
2784  {
2785  om_->stream(Debug) << "\nSeeking a shift for theta[0]=" << thetaMag[0] << std::endl;
2786 
2787  // This is my adjustment. If all eigenvalues are in a single cluster, it's probably a bad idea to shift the smallest one.
2788  // If all of your eigenvalues are in one cluster, it's either way to early to shift or your subspace is too small
2789  if((theta_[0] > 0 && posSafeToShift_) || (theta_[0] < 0 && negSafeToShift_) || considerClusters_ == false)
2790  {
2791  // Initialize with a conservative shift, either the biggest safe shift or the eigenvalue adjusted by its cluster's residual
2792  ritzShifts_[0] = std::max(largestSafeShift_,thetaMag[0]-clusterResids[0]);
2793 
2794  om_->stream(Debug) << "Initializing with a conservative shift, either the most positive converged eigenvalue ("
2795  << largestSafeShift_ << ") or the eigenvalue adjusted by the residual (" << thetaMag[0] << "-"
2796  << clusterResids[0] << ").\n";
2797 
2798  // If this eigenvalue is NOT in a cluster, do an aggressive shift
2799  if(R2norms_[0] == clusterResids[0])
2800  {
2801  ritzShifts_[0] = thetaMag[0];
2802  om_->stream(Debug) << "Since this eigenvalue is NOT in a cluster, we can use the eigenvalue itself as a shift: ritzShifts[0]=" << ritzShifts_[0] << std::endl;
2803  }
2804  else
2805  om_->stream(Debug) << "This eigenvalue is in a cluster, so it would not be safe to use the eigenvalue itself as a shift\n";
2806  }
2807  else
2808  {
2809  if(largestSafeShift_ > std::abs(ritzShifts_[0]))
2810  {
2811  om_->stream(Debug) << "Initializing with a conservative shift...the most positive converged eigenvalue: " << largestSafeShift_ << std::endl;
2812  ritzShifts_[0] = largestSafeShift_;
2813  }
2814  else
2815  om_->stream(Debug) << "Using the previous value of ritzShifts[0]=" << ritzShifts_[0];
2816 
2817  }
2818 
2819  om_->stream(Debug) << "ritzShifts[0]=" << ritzShifts_[0] << std::endl;
2820 
2821  if(useMultipleShifts_)
2822  {
2824  // Compute shifts for other eigenvalues
2825  for(int i=1; i < blockSize_; i++)
2826  {
2827  om_->stream(Debug) << "\nSeeking a shift for theta[" << i << "]=" << thetaMag[i] << std::endl;
2828 
2829  // If the previous shift was aggressive and we are not in a cluster, do an aggressive shift
2830  if(ritzShifts_[i-1] == thetaMag[i-1] && i < blockSize_-1 && thetaMag[i] < thetaMag[i+1] - clusterResids[i+1])
2831  {
2832  ritzShifts_[i] = thetaMag[i];
2833  om_->stream(Debug) << "Using an aggressive shift: ritzShifts_[" << i << "]=" << ritzShifts_[i] << std::endl;
2834  }
2835  else
2836  {
2837  if(ritzShifts_[0] > std::abs(ritzShifts_[i]))
2838  {
2839  om_->stream(Debug) << "It was unsafe to use the aggressive shift. Choose the shift used by theta[0]="
2840  << thetaMag[0] << ": ritzShifts[0]=" << ritzShifts_[0] << std::endl;
2841 
2842  // Choose a conservative shift, that of the smallest positive eigenvalue
2843  ritzShifts_[i] = ritzShifts_[0];
2844  }
2845  else
2846  om_->stream(Debug) << "It was unsafe to use the aggressive shift. We will use the shift from the previous iteration: " << ritzShifts_[i] << std::endl;
2847 
2848  om_->stream(Debug) << "Check whether any less conservative shifts would work (such as the biggest eigenvalue outside of the cluster, namely theta[ell] < "
2849  << thetaMag[i] << "-" << clusterResids[i] << " (" << thetaMag[i] - clusterResids[i] << ")\n";
2850 
2851  // If possible, choose a less conservative shift, that of the biggest eigenvalue outside of the cluster
2852  for(int ell=0; ell < i; ell++)
2853  {
2854  if(thetaMag[ell] < thetaMag[i] - clusterResids[i])
2855  {
2856  ritzShifts_[i] = thetaMag[ell];
2857  om_->stream(Debug) << "ritzShifts_[" << i << "]=" << ritzShifts_[i] << " is valid\n";
2858  }
2859  else
2860  break;
2861  }
2862  } // end else
2863 
2864  om_->stream(Debug) << "ritzShifts[" << i << "]=" << ritzShifts_[i] << std::endl;
2865  } // end for
2866  } // end if(useMultipleShifts_)
2867  else
2868  {
2869  for(int i=1; i<blockSize_; i++)
2870  ritzShifts_[i] = ritzShifts_[0];
2871  }
2872  } // end if(howToShift_ == "Adjusted Ritz Values")
2873  } // end if(whenToShift_ == "Always" || (whenToShift_ == "After Trace Levels" && traceHasLeveled))
2874 
2875  // Set the correct sign
2876  for(int i=0; i<blockSize_; i++)
2877  {
2878  if(theta_[i] < 0)
2879  ritzShifts_[i] = -abs(ritzShifts_[i]);
2880  else
2881  ritzShifts_[i] = abs(ritzShifts_[i]);
2882  }
2883  }
2884 
2886  template <class ScalarType, class MV, class OP>
2887  std::vector<ScalarType> TraceMinBase<ScalarType,MV,OP>::computeTol()
2888  {
2889  ScalarType temp1;
2890  std::vector<ScalarType> tolerances(blockSize_);
2891 
2892  for(int i=0; i < blockSize_-1; i++)
2893  {
2894  if(std::abs(theta_[blockSize_-1]) != std::abs(ritzShifts_[i]))
2895  temp1 = std::abs(theta_[i]-ritzShifts_[i])/std::abs(std::abs(theta_[blockSize_-1])-std::abs(ritzShifts_[i]));
2896  else
2897  temp1 = ZERO;
2898 
2899  // TODO: The min and max tolerances should not be hard coded
2900  // Neither should the maximum number of iterations
2901  tolerances[i] = std::min(temp1*temp1,0.5);
2902  }
2903 
2904  if(blockSize_ > 1)
2905  tolerances[blockSize_-1] = tolerances[blockSize_-2];
2906 
2907  return tolerances;
2908  }
2909 
2911  template <class ScalarType, class MV, class OP>
2912  void TraceMinBase<ScalarType,MV,OP>::solveSaddlePointProblem(RCP<MV> Delta)
2913  {
2914 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
2915  Teuchos::TimeMonitor lcltimer( *timerSaddle_ );
2916 #endif
2917 
2918  // This case can arise when looking for the largest eigenpairs
2919  if(Op_ == Teuchos::null)
2920  {
2921  // dense solver
2923 
2924  // Schur complement
2926 
2927  if(computeAllRes_)
2928  {
2929  // Get the valid indices of X
2930  std::vector<int> curind(blockSize_);
2931  for(int i=0; i<blockSize_; i++)
2932  curind[i] = i;
2933 
2934  // Get a view of MX
2935  RCP<const MV> lclMX = MVT::CloneView(*MX_, curind);
2936 
2937  // form S = X' M^2 X
2938  MVT::MvTransMv(ONE,*lclMX,*lclMX,*lclS);
2939 
2940  // compute the inverse of S
2941  My_Solver.setMatrix(lclS);
2942  My_Solver.invert();
2943 
2944  // Delta = X - MX/S
2945  RCP<const MV> lclX = MVT::CloneView(*X_, curind);
2946  MVT::Assign(*lclX,*Delta);
2947  MVT::MvTimesMatAddMv( -ONE, *lclMX, *lclS, ONE, *Delta);
2948  }
2949  else
2950  {
2951  // form S = X' M^2 X
2952  MVT::MvTransMv(ONE,*MX_,*MX_,*lclS);
2953 
2954  // compute the inverse of S
2955  My_Solver.setMatrix(lclS);
2956  My_Solver.invert();
2957 
2958  // Delta = X - MX/S
2959  MVT::Assign(*X_,*Delta);
2960  MVT::MvTimesMatAddMv( -ONE, *MX_, *lclS, ONE, *Delta);
2961  }
2962  }
2963  else
2964  {
2965  std::vector<int> order(curDim_);
2966  std::vector<ScalarType> tempvec(blockSize_);
2967 // RCP<BasicSort<MagnitudeType> > sorter = rcp( new BasicSort<MagnitudeType>("SR") );
2968 
2969  // Stores the residual of each CLUSTER of eigenvalues
2970  std::vector<ScalarType> clusterResids;
2971 
2972 /* // Sort the eigenvalues in ascending order for the Ritz shift selection
2973  sorter->sort(theta_, Teuchos::rcpFromRef(order), curDim_); // don't catch exception
2974 
2975  // Apply the same ordering to the residual norms
2976  getRes2Norms();
2977  for (int i=0; i<blockSize_; i++)
2978  tempvec[i] = R2norms_[order[i]];
2979  R2norms_ = tempvec;*/
2980 
2981  // Compute the residual of each CLUSTER of eigenvalues
2982  // This is important for selecting the Ritz shifts
2983  clusterResids = getClusterResids();
2984 
2985 /* // Sort the eigenvalues based on what the user wanted
2986  sm_->sort(theta_, Teuchos::rcpFromRef(order), blockSize_);
2987 
2988  // Apply the same ordering to the residual norms and cluster residuals
2989  for (int i=0; i<blockSize_; i++)
2990  tempvec[i] = R2norms_[order[i]];
2991  R2norms_ = tempvec;
2992 
2993  for (int i=0; i<blockSize_; i++)
2994  tempvec[i] = clusterResids[order[i]];
2995  clusterResids = tempvec;*/
2996 
2997  // Compute the Ritz shifts
2998  computeRitzShifts(clusterResids);
2999 
3000  // Compute the tolerances for the inner solves
3001  std::vector<ScalarType> tolerances = computeTol();
3002 
3003  for(int i=0; i<blockSize_; i++)
3004  {
3005  om_->stream(IterationDetails) << "Choosing Ritz shifts...theta[" << i << "]="
3006  << theta_[i] << ", resids[" << i << "]=" << R2norms_[i] << ", clusterResids[" << i << "]=" << clusterResids[i]
3007  << ", ritzShifts[" << i << "]=" << ritzShifts_[i] << ", and tol[" << i << "]=" << tolerances[i] << std::endl;
3008  }
3009 
3010  // Set the Ritz shifts for the solver
3011  ritzOp_->setRitzShifts(ritzShifts_);
3012 
3013  // Set the inner stopping tolerance
3014  // This uses the Ritz values to determine when to stop
3015  ritzOp_->setInnerTol(tolerances);
3016 
3017  // Solve the saddle point problem
3018  if(saddleSolType_ == PROJECTED_KRYLOV_SOLVER)
3019  {
3020  if(Prec_ != Teuchos::null)
3021  solveSaddleProjPrec(Delta);
3022  else
3023  solveSaddleProj(Delta);
3024  }
3025  else if(saddleSolType_ == SCHUR_COMPLEMENT_SOLVER)
3026  {
3027  if(Z_ == Teuchos::null || MVT::GetNumberVecs(*Z_) != blockSize_)
3028  {
3029  // We do NOT want Z to be 0, because that could result in stagnation
3030  // I know it's tempting to take out the MvRandom, but seriously, don't do it.
3031  Z_ = MVT::Clone(*X_,blockSize_);
3032  MVT::MvRandom(*Z_);
3033  }
3034  solveSaddleSchur(Delta);
3035  }
3036  else if(saddleSolType_ == BD_PREC_MINRES)
3037  {
3038  solveSaddleBDPrec(Delta);
3039 // Delta->describe(*(Teuchos::VerboseObjectBase::getDefaultOStream()),Teuchos::VERB_EXTREME);
3040  }
3041  else if(saddleSolType_ == HSS_PREC_GMRES)
3042  {
3043  solveSaddleHSSPrec(Delta);
3044  }
3045  else
3046  std::cout << "Invalid saddle solver type\n";
3047  }
3048  }
3049 
3051  // Solve the saddle point problem using projected minres
3052  // TODO: We should be able to choose KX or -R as RHS.
3053  template <class ScalarType, class MV, class OP>
3054  void TraceMinBase<ScalarType,MV,OP>::solveSaddleProj (RCP<MV> Delta) const
3055  {
3057 
3058  if(computeAllRes_)
3059  {
3060  // Get the valid indices of X
3061  std::vector<int> curind(blockSize_);
3062  for(int i=0; i<blockSize_; i++)
3063  curind[i] = i;
3064 
3065  RCP<const MV> locMX = MVT::CloneView(*MX_, curind);
3066 
3067  // We should really project out the converged eigenvectors too
3068  if(projectAllVecs_)
3069  {
3070  if(projectLockedVecs_ && numAuxVecs_ > 0)
3071  projOp = rcp( new TraceMinProjRitzOp<ScalarType,MV,OP>(ritzOp_,locMX,orthman_,auxVecs_) );
3072  else
3073  projOp = rcp( new TraceMinProjRitzOp<ScalarType,MV,OP>(ritzOp_,locMX,orthman_) );
3074  }
3075  else
3076  projOp = rcp( new TraceMinProjRitzOp<ScalarType,MV,OP>(ritzOp_,locMX) );
3077 
3078  // Remember, Delta0 must equal 0
3079  // This ensures B-orthogonality between Delta and X
3080  MVT::MvInit(*Delta);
3081 
3082  if(useRHSR_)
3083  {
3084  RCP<const MV> locR = MVT::CloneView(*R_, curind);
3085  projOp->ApplyInverse(*locR, *Delta);
3086  }
3087  else
3088  {
3089  RCP<const MV> locKX = MVT::CloneView(*KX_, curind);
3090  projOp->ApplyInverse(*locKX, *Delta);
3091  }
3092  }
3093  else
3094  {
3095  // We should really project out the converged eigenvectors too
3096  if(projectAllVecs_)
3097  {
3098  if(projectLockedVecs_ && numAuxVecs_ > 0)
3099  projOp = rcp( new TraceMinProjRitzOp<ScalarType,MV,OP>(ritzOp_,MX_,orthman_,auxVecs_) );
3100  else
3101  projOp = rcp( new TraceMinProjRitzOp<ScalarType,MV,OP>(ritzOp_,MX_,orthman_) );
3102  }
3103  else
3104  projOp = rcp( new TraceMinProjRitzOp<ScalarType,MV,OP>(ritzOp_,MX_) );
3105 
3106  // Remember, Delta0 must equal 0
3107  // This ensures B-orthogonality between Delta and X
3108  MVT::MvInit(*Delta);
3109 
3110  if(useRHSR_) {
3111  projOp->ApplyInverse(*R_, *Delta);
3112  }
3113  else {
3114  projOp->ApplyInverse(*KX_, *Delta);
3115  }
3116  }
3117  }
3118 
3120  // TODO: Fix preconditioning
3121  template <class ScalarType, class MV, class OP>
3122  void TraceMinBase<ScalarType,MV,OP>::solveSaddleProjPrec (RCP<MV> Delta) const
3123  {
3124  // If we don't have Belos installed, we can't use TraceMinProjRitzOpWithPrec
3125  // Of course, this problem will be detected in the constructor and an exception will be thrown
3126  // This is only here to make sure the code compiles properly
3127 #ifdef HAVE_ANASAZI_BELOS
3129 
3130  if(computeAllRes_)
3131  {
3132  int dimension;
3133  if(projectAllVecs_)
3134  dimension = curDim_;
3135  else
3136  dimension = blockSize_;
3137 
3138  // Get the valid indices of X
3139  std::vector<int> curind(dimension);
3140  for(int i=0; i<dimension; i++)
3141  curind[i] = i;
3142 
3143  RCP<const MV> locMX = MVT::CloneView(*MX_, curind);
3144 
3145  // We should really project out the converged eigenvectors too
3146  if(projectAllVecs_)
3147  {
3148  if(projectLockedVecs_ && numAuxVecs_ > 0)
3149  projOp = rcp( new TraceMinProjRitzOpWithPrec<ScalarType,MV,OP>(ritzOp_,locMX,orthman_,auxVecs_) );
3150  else
3151  projOp = rcp( new TraceMinProjRitzOpWithPrec<ScalarType,MV,OP>(ritzOp_,locMX,orthman_) );
3152  }
3153  else
3154  projOp = rcp( new TraceMinProjRitzOpWithPrec<ScalarType,MV,OP>(ritzOp_,locMX) );
3155 
3156  // Remember, Delta0 must equal 0
3157  // This ensures B-orthogonality between Delta and X
3158  MVT::MvInit(*Delta);
3159 
3160  std::vector<int> dimind(blockSize_);
3161  for(int i=0; i<blockSize_; i++)
3162  dimind[i] = i;
3163 
3164  if(useRHSR_)
3165  {
3166  RCP<const MV> locR = MVT::CloneView(*R_, dimind);
3167  projOp->ApplyInverse(*locR, *Delta);
3168  MVT::MvScale(*Delta, -ONE);
3169  }
3170  else
3171  {
3172  RCP<const MV> locKX = MVT::CloneView(*KX_, dimind);
3173  projOp->ApplyInverse(*locKX, *Delta);
3174  }
3175  }
3176  else
3177  {
3178  // We should really project out the converged eigenvectors too
3179  if(projectAllVecs_)
3180  {
3181  if(projectLockedVecs_ && numAuxVecs_ > 0)
3182  projOp = rcp( new TraceMinProjRitzOpWithPrec<ScalarType,MV,OP>(ritzOp_,MX_,orthman_,auxVecs_) );
3183  else
3184  projOp = rcp( new TraceMinProjRitzOpWithPrec<ScalarType,MV,OP>(ritzOp_,MX_,orthman_) );
3185  }
3186  else
3187  projOp = rcp( new TraceMinProjRitzOpWithPrec<ScalarType,MV,OP>(ritzOp_,MX_) );
3188 
3189  // Remember, Delta0 must equal 0
3190  // This ensures B-orthogonality between Delta and X
3191  MVT::MvInit(*Delta);
3192 
3193  if(useRHSR_)
3194  {
3195  projOp->ApplyInverse(*R_, *Delta);
3196  MVT::MvScale(*Delta,-ONE);
3197  }
3198  else
3199  projOp->ApplyInverse(*KX_, *Delta);
3200  }
3201 #endif
3202  }
3203 
3205  // TODO: We can hold the Schur complement constant in later iterations
3206  // TODO: Make sure we're using the preconditioner correctly
3207  template <class ScalarType, class MV, class OP>
3208  void TraceMinBase<ScalarType,MV,OP>::solveSaddleSchur (RCP<MV> Delta) const
3209  {
3210  // dense solver
3212 
3215 
3216  if(computeAllRes_)
3217  {
3218  // Get the valid indices of X
3219  std::vector<int> curind(blockSize_);
3220  for(int i=0; i<blockSize_; i++)
3221  curind[i] = i;
3222 
3223  // Z = K \ MX
3224  // Why would I have wanted to set Z <- X? I want to leave Z's previous value alone...
3225  RCP<const MV> lclMX = MVT::CloneView(*MX_, curind);
3226 
3227 #ifdef USE_APPLY_INVERSE
3228  Op_->ApplyInverse(*lclMX,*Z_);
3229 #else
3230  ritzOp_->ApplyInverse(*lclMX,*Z_);
3231 #endif
3232 
3233  // form S = X' M Z
3234  MVT::MvTransMv(ONE,*Z_,*lclMX,*lclS);
3235 
3236  // solve S L = I
3237  My_Solver.setMatrix(lclS);
3238  My_Solver.invert();
3239  lclL = lclS;
3240 
3241  // Delta = X - Z L
3242  RCP<const MV> lclX = MVT::CloneView(*X_, curind);
3243  MVT::Assign(*lclX,*Delta);
3244  MVT::MvTimesMatAddMv( -ONE, *Z_, *lclL, ONE, *Delta);
3245  }
3246  else
3247  {
3248  // Z = K \ MX
3249 #ifdef USE_APPLY_INVERSE
3250  Op_->ApplyInverse(*MX_,*Z_);
3251 #else
3252  ritzOp_->ApplyInverse(*MX_,*Z_);
3253 #endif
3254 
3255  // form S = X' M Z
3256  MVT::MvTransMv(ONE,*Z_,*MX_,*lclS);
3257 
3258  // solve S L = I
3259  My_Solver.setMatrix(lclS);
3260  My_Solver.invert();
3261  lclL = lclS;
3262 
3263  // Delta = X - Z L
3264  MVT::Assign(*X_,*Delta);
3265  MVT::MvTimesMatAddMv( -ONE, *Z_, *lclL, ONE, *Delta);
3266  }
3267  }
3268 
3269 
3271  // TODO: We can hold the Schur complement constant in later iterations
3272  template <class ScalarType, class MV, class OP>
3273  void TraceMinBase<ScalarType,MV,OP>::solveSaddleBDPrec (RCP<MV> Delta) const
3274  {
3275  RCP<MV> locKX, locMX;
3276  if(computeAllRes_)
3277  {
3278  std::vector<int> curind(blockSize_);
3279  for(int i=0; i<blockSize_; i++)
3280  curind[i] = i;
3281  locKX = MVT::CloneViewNonConst(*KX_, curind);
3282  locMX = MVT::CloneViewNonConst(*MX_, curind);
3283  }
3284  else
3285  {
3286  locKX = KX_;
3287  locMX = MX_;
3288  }
3289 
3290  // Create the operator [A BX; X'B 0]
3291  RCP<saddle_op_type> sadOp = rcp(new saddle_op_type(ritzOp_,locMX));
3292 
3293  // Create the RHS [AX; 0]
3294  RCP<saddle_container_type> sadRHS = rcp(new saddle_container_type(locKX));
3295 
3296 // locKX->describe(*(Teuchos::VerboseObjectBase::getDefaultOStream()),Teuchos::VERB_EXTREME);
3297 
3298 // locMX->describe(*(Teuchos::VerboseObjectBase::getDefaultOStream()),Teuchos::VERB_EXTREME);
3299 
3300  // Create the solution vector [Delta; L]
3301  MVT::MvInit(*Delta);
3302  RCP<saddle_container_type> sadSol = rcp(new saddle_container_type(Delta));
3303 
3304  // Create a minres solver
3306  if(Prec_ != Teuchos::null)
3307  {
3308  RCP<saddle_op_type> sadPrec = rcp(new saddle_op_type(ritzOp_->getPrec(),locMX,BD_PREC));
3309  sadSolver = rcp(new PseudoBlockMinres<ScalarType,saddle_container_type,saddle_op_type>(sadOp, sadPrec));
3310  }
3311  else {
3312  sadSolver = rcp(new PseudoBlockMinres<ScalarType,saddle_container_type,saddle_op_type>(sadOp));
3313  }
3314 
3315  // Set the tolerance for the minres solver
3316  std::vector<ScalarType> tol;
3317  ritzOp_->getInnerTol(tol);
3318  sadSolver->setTol(tol);
3319 
3320  // Set the maximum number of iterations
3321  sadSolver->setMaxIter(ritzOp_->getMaxIts());
3322 
3323  // Set the solution vector
3324  sadSolver->setSol(sadSol);
3325 
3326  // Set the RHS
3327  sadSolver->setRHS(sadRHS);
3328 
3329  // Solve the saddle point problem
3330  sadSolver->solve();
3331  }
3332 
3333 
3334 
3336  // TODO: We can hold the Schur complement constant in later iterations
3337  template <class ScalarType, class MV, class OP>
3338  void TraceMinBase<ScalarType,MV,OP>::solveSaddleHSSPrec (RCP<MV> Delta) const
3339  {
3340 #ifdef HAVE_ANASAZI_BELOS
3341  typedef Belos::LinearProblem<ScalarType,saddle_container_type,saddle_op_type> LP;
3342  typedef Belos::PseudoBlockGmresSolMgr<ScalarType,saddle_container_type,saddle_op_type> GmresSolMgr;
3343 
3344  RCP<MV> locKX, locMX;
3345  if(computeAllRes_)
3346  {
3347  std::vector<int> curind(blockSize_);
3348  for(int i=0; i<blockSize_; i++)
3349  curind[i] = i;
3350  locKX = MVT::CloneViewNonConst(*KX_, curind);
3351  locMX = MVT::CloneViewNonConst(*MX_, curind);
3352  }
3353  else
3354  {
3355  locKX = KX_;
3356  locMX = MX_;
3357  }
3358 
3359  // Create the operator [A BX; X'B 0]
3360  RCP<saddle_op_type> sadOp = rcp(new saddle_op_type(ritzOp_,locMX,NONSYM));
3361 
3362  // Create the RHS [AX; 0]
3363  RCP<saddle_container_type> sadRHS = rcp(new saddle_container_type(locKX));
3364 
3365  // Create the solution vector [Delta; L]
3366  MVT::MvInit(*Delta);
3367  RCP<saddle_container_type> sadSol = rcp(new saddle_container_type(Delta));
3368 
3369  // Create a parameter list for the gmres solver
3371 
3372  // Set the tolerance for the gmres solver
3373  std::vector<ScalarType> tol;
3374  ritzOp_->getInnerTol(tol);
3375  pl->set("Convergence Tolerance",tol[0]);
3376 
3377  // Set the maximum number of iterations
3378  // TODO: Come back to this
3379  pl->set("Maximum Iterations", ritzOp_->getMaxIts());
3380  pl->set("Num Blocks", ritzOp_->getMaxIts());
3381 
3382  // Set the block size
3383  // TODO: Come back to this
3384  // TODO: This breaks the code right now, presumably because of a MVT cloneview issue.
3385  pl->set("Block Size", blockSize_);
3386 
3387  // Set the verbosity of gmres
3388 // pl->set("Verbosity", Belos::IterationDetails + Belos::StatusTestDetails + Belos::Debug);
3389 // pl->set("Output Frequency", 1);
3390 
3391  // Create the linear problem
3392  RCP<LP> problem = rcp(new LP(sadOp,sadSol,sadRHS));
3393 
3394  // Set the preconditioner
3395  if(Prec_ != Teuchos::null)
3396  {
3397  RCP<saddle_op_type> sadPrec = rcp(new saddle_op_type(ritzOp_->getPrec(),locMX,HSS_PREC,alpha_));
3398  problem->setLeftPrec(sadPrec);
3399  }
3400 
3401  // Set the problem
3402  problem->setProblem();
3403 
3404  // Create a minres solver
3405  RCP<GmresSolMgr> sadSolver = rcp(new GmresSolMgr(problem,pl)) ;
3406 
3407  // Solve the saddle point problem
3408  sadSolver->solve();
3409 #else
3410  std::cout << "No Belos. This is bad\n";
3411 #endif
3412  }
3413 
3414 
3415 
3417  // Compute KK := V'KV
3418  // We only compute the NEW elements
3419  template <class ScalarType, class MV, class OP>
3420  void TraceMinBase<ScalarType,MV,OP>::computeKK()
3421  {
3422  // Get the valid indices of V
3423  std::vector<int> curind(curDim_);
3424  for(int i=0; i<curDim_; i++)
3425  curind[i] = i;
3426 
3427  // Get a pointer to the valid parts of V
3428  RCP<const MV> lclV = MVT::CloneView(*V_,curind);
3429 
3430  // Get the valid indices of KV
3431  curind.resize(blockSize_);
3432  for(int i=0; i<blockSize_; i++)
3433  curind[i] = curDim_-blockSize_+i;
3434  RCP<const MV> lclKV = MVT::CloneView(*KV_,curind);
3435 
3436  // Get a pointer to the valid part of KK
3438  rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,blockSize_,0,curDim_-blockSize_) );
3439 
3440  // KK := V'KV
3441  MVT::MvTransMv(ONE,*lclV,*lclKV,*lclKK);
3442 
3443  // We only constructed the upper triangular part of the matrix, but that's okay because KK is symmetric!
3444  for(int r=0; r<curDim_; r++)
3445  {
3446  for(int c=0; c<r; c++)
3447  {
3448  (*KK_)(r,c) = (*KK_)(c,r);
3449  }
3450  }
3451  }
3452 
3454  // Compute the eigenpairs of KK, i.e. the Ritz pairs
3455  template <class ScalarType, class MV, class OP>
3456  void TraceMinBase<ScalarType,MV,OP>::computeRitzPairs()
3457  {
3458  // Get a pointer to the valid part of KK
3460  rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,curDim_) );
3461 
3462  // Get a pointer to the valid part of ritzVecs
3464  rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,curDim_) );
3465 
3466  // Compute Ritz pairs from KK
3467  {
3468 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
3469  Teuchos::TimeMonitor lcltimer( *timerDS_ );
3470 #endif
3471  int rank = curDim_;
3472  Utils::directSolver(curDim_, *lclKK, Teuchos::null, *lclRV, theta_, rank, 10);
3473  // we want all ritz values back
3474  // TODO: This probably should not be an ortho failure
3475  TEUCHOS_TEST_FOR_EXCEPTION(rank != curDim_,TraceMinBaseOrthoFailure,
3476  "Anasazi::TraceMinBase::computeRitzPairs(): Failed to compute all eigenpairs of KK.");
3477  }
3478 
3479  // Sort ritz pairs
3480  {
3481 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
3482  Teuchos::TimeMonitor lcltimer( *timerSortEval_ );
3483 #endif
3484 
3485  std::vector<int> order(curDim_);
3486  //
3487  // sort the first curDim_ values in theta_
3488  if(useHarmonic_)
3489  {
3491  sm.sort(theta_, Teuchos::rcpFromRef(order), curDim_);
3492  }
3493  else
3494  {
3495  sm_->sort(theta_, Teuchos::rcpFromRef(order), curDim_); // don't catch exception
3496  }
3497  //
3498  // apply the same ordering to the primitive ritz vectors
3499  Utils::permuteVectors(order,*lclRV);
3500  }
3501  }
3502 
3504  // Compute X := V evecs
3505  template <class ScalarType, class MV, class OP>
3506  void TraceMinBase<ScalarType,MV,OP>::computeX()
3507  {
3508 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
3509  Teuchos::TimeMonitor lcltimer( *timerLocal_ );
3510 #endif
3511 
3512  // Get the valid indices of V
3513  std::vector<int> curind(curDim_);
3514  for(int i=0; i<curDim_; i++)
3515  curind[i] = i;
3516 
3517  // Get a pointer to the valid parts of V
3518  RCP<const MV> lclV = MVT::CloneView(*V_,curind);
3519 
3520  if(computeAllRes_)
3521  {
3522  // Capture the relevant eigenvectors
3524  rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,curDim_) );
3525 
3526  // X <- lclV*S
3527  RCP<MV> lclX = MVT::CloneViewNonConst(*X_,curind);
3528  MVT::MvTimesMatAddMv( ONE, *lclV, *relevantEvecs, ZERO, *lclX );
3529  }
3530  else
3531  {
3532  // Capture the relevant eigenvectors
3534  rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,blockSize_) );
3535 
3536  // X <- lclV*S
3537  MVT::MvTimesMatAddMv( ONE, *lclV, *relevantEvecs, ZERO, *X_ );
3538  }
3539  }
3540 
3542  // Compute KX := KV evecs and (if necessary) MX := MV evecs
3543  template <class ScalarType, class MV, class OP>
3544  void TraceMinBase<ScalarType,MV,OP>::updateKXMX()
3545  {
3546 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
3547  Teuchos::TimeMonitor lcltimer( *timerLocal_ );
3548 #endif
3549 
3550  // Get the valid indices of V
3551  std::vector<int> curind(curDim_);
3552  for(int i=0; i<curDim_; i++)
3553  curind[i] = i;
3554 
3555  // Get pointers to the valid parts of V, KV, and MV (if necessary)
3556  RCP<const MV> lclV = MVT::CloneView(*V_,curind);
3557  RCP<const MV> lclKV = MVT::CloneView(*KV_,curind);
3558 
3559  if(computeAllRes_)
3560  {
3561  // Capture the relevant eigenvectors
3563  rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,curDim_) );
3564 
3565  // Update KX and MX
3566  RCP<MV> lclKX = MVT::CloneViewNonConst(*KX_,curind);
3567  MVT::MvTimesMatAddMv( ONE, *lclKV, *relevantEvecs, ZERO, *lclKX );
3568  if(hasM_)
3569  {
3570  RCP<const MV> lclMV = MVT::CloneView(*MV_,curind);
3571  RCP<MV> lclMX = MVT::CloneViewNonConst(*MX_,curind);
3572  MVT::MvTimesMatAddMv( ONE, *lclMV, *relevantEvecs, ZERO, *lclMX );
3573  }
3574  }
3575  else
3576  {
3577  // Capture the relevant eigenvectors
3579  rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*ritzVecs_,curDim_,blockSize_) );
3580 
3581  // Update KX and MX
3582  MVT::MvTimesMatAddMv( ONE, *lclKV, *relevantEvecs, ZERO, *KX_ );
3583  if(hasM_)
3584  {
3585  RCP<const MV> lclMV = MVT::CloneView(*MV_,curind);
3586  MVT::MvTimesMatAddMv( ONE, *lclMV, *relevantEvecs, ZERO, *MX_ );
3587  }
3588  }
3589  }
3590 
3592  // Update the residual R := KX - MX*T
3593  template <class ScalarType, class MV, class OP>
3594  void TraceMinBase<ScalarType,MV,OP>::updateResidual () {
3595 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
3596  Teuchos::TimeMonitor lcltimer( *timerCompRes_ );
3597 #endif
3598 
3599  if(computeAllRes_)
3600  {
3601  // Get the valid indices of X
3602  std::vector<int> curind(curDim_);
3603  for(int i=0; i<curDim_; i++)
3604  curind[i] = i;
3605 
3606  // Holds MX*T
3607  RCP<MV> MXT = MVT::CloneCopy(*MX_,curind);
3608 
3609  // Holds the relevant part of theta
3610  std::vector<ScalarType> locTheta(curDim_);
3611  for(int i=0; i<curDim_; i++)
3612  locTheta[i] = theta_[i];
3613 
3614  // Compute MX*T
3615  MVT::MvScale(*MXT,locTheta);
3616 
3617  // form R <- KX - MX*T
3618  RCP<const MV> locKX = MVT::CloneView(*KX_,curind);
3619  RCP<MV> locR = MVT::CloneViewNonConst(*R_,curind);
3620  MVT::MvAddMv(ONE,*locKX,-ONE,*MXT,*locR);
3621  }
3622  else
3623  {
3624  // Holds MX*T
3625  RCP<MV> MXT = MVT::CloneCopy(*MX_);
3626 
3627  // Holds the relevant part of theta
3628  std::vector<ScalarType> locTheta(blockSize_);
3629  for(int i=0; i<blockSize_; i++)
3630  locTheta[i] = theta_[i];
3631 
3632  // Compute MX*T
3633  MVT::MvScale(*MXT,locTheta);
3634 
3635  // form R <- KX - MX*T
3636  MVT::MvAddMv(ONE,*KX_,-ONE,*MXT,*R_);
3637  }
3638 
3639  // R has been updated; mark the norms as out-of-date
3640  Rnorms_current_ = false;
3641  R2norms_current_ = false;
3642  }
3643 
3644 
3646  // Check accuracy, orthogonality, and other debugging stuff
3647  //
3648  // bools specify which tests we want to run (instead of running more than we actually care about)
3649  //
3650  // we don't bother checking the following because they are computed explicitly:
3651  // H == Prec*R
3652  // KH == K*H
3653  //
3654  //
3655  // checkV : V orthonormal
3656  // orthogonal to auxvecs
3657  // checkX : X orthonormal
3658  // orthogonal to auxvecs
3659  // checkMX: check MX == M*X
3660  // checkKX: check KX == K*X
3661  // checkH : H orthonormal
3662  // orthogonal to V and H and auxvecs
3663  // checkMH: check MH == M*H
3664  // checkR : check R orthogonal to X
3665  // checkQ : check that auxiliary vectors are actually orthonormal
3666  // checkKK: check that KK is symmetric in memory
3667  //
3668  // TODO:
3669  // add checkTheta
3670  //
3671  template <class ScalarType, class MV, class OP>
3672  std::string TraceMinBase<ScalarType,MV,OP>::accuracyCheck( const CheckList &chk, const std::string &where ) const
3673  {
3674  using std::endl;
3675 
3676  std::stringstream os;
3677  os.precision(2);
3678  os.setf(std::ios::scientific, std::ios::floatfield);
3679 
3680  os << " Debugging checks: iteration " << iter_ << where << endl;
3681 
3682  // V and friends
3683  std::vector<int> lclind(curDim_);
3684  for (int i=0; i<curDim_; ++i) lclind[i] = i;
3685  RCP<const MV> lclV;
3686  if (initialized_) {
3687  lclV = MVT::CloneView(*V_,lclind);
3688  }
3689  if (chk.checkV && initialized_) {
3690  MagnitudeType err = orthman_->orthonormError(*lclV);
3691  os << " >> Error in V^H M V == I : " << err << endl;
3692  for (Array_size_type i=0; i<auxVecs_.size(); ++i) {
3693  err = orthman_->orthogError(*lclV,*auxVecs_[i]);
3694  os << " >> Error in V^H M Q[" << i << "] == 0 : " << err << endl;
3695  }
3696  }
3697 
3698  // X and friends
3699  RCP<const MV> lclX;
3700  if(initialized_)
3701  {
3702  if(computeAllRes_)
3703  lclX = MVT::CloneView(*X_,lclind);
3704  else
3705  lclX = X_;
3706  }
3707 
3708  if (chk.checkX && initialized_) {
3709  MagnitudeType err = orthman_->orthonormError(*lclX);
3710  os << " >> Error in X^H M X == I : " << err << endl;
3711  for (Array_size_type i=0; i<auxVecs_.size(); ++i) {
3712  err = orthman_->orthogError(*lclX,*auxVecs_[i]);
3713  os << " >> Error in X^H M Q[" << i << "] == 0 : " << err << endl;
3714  }
3715  }
3716  if (chk.checkMX && hasM_ && initialized_) {
3717  RCP<const MV> lclMX;
3718  if(computeAllRes_)
3719  lclMX = MVT::CloneView(*MX_,lclind);
3720  else
3721  lclMX = MX_;
3722 
3723  MagnitudeType err = Utils::errorEquality(*lclX, *lclMX, MOp_);
3724  os << " >> Error in MX == M*X : " << err << endl;
3725  }
3726  if (Op_ != Teuchos::null && chk.checkKX && initialized_) {
3727  RCP<const MV> lclKX;
3728  if(computeAllRes_)
3729  lclKX = MVT::CloneView(*KX_,lclind);
3730  else
3731  lclKX = KX_;
3732 
3733  MagnitudeType err = Utils::errorEquality(*lclX, *lclKX, Op_);
3734  os << " >> Error in KX == K*X : " << err << endl;
3735  }
3736 
3737  // KK
3738  if (chk.checkKK && initialized_) {
3739  Teuchos::SerialDenseMatrix<int,ScalarType> curKK(curDim_,curDim_);
3740  if(Op_ != Teuchos::null) {
3741  RCP<MV> lclKV = MVT::Clone(*V_,curDim_);
3742  OPT::Apply(*Op_,*lclV,*lclKV);
3743  MVT::MvTransMv(ONE,*lclV,*lclKV,curKK);
3744  }
3745  else {
3746  MVT::MvTransMv(ONE,*lclV,*lclV,curKK);
3747  }
3748  Teuchos::SerialDenseMatrix<int,ScalarType> subKK(Teuchos::View,*KK_,curDim_,curDim_);
3749  curKK -= subKK;
3750  os << " >> Error in V^H K V == KK : " << curKK.normFrobenius() << endl;
3751 
3752  Teuchos::SerialDenseMatrix<int,ScalarType> SDMerr(curDim_,curDim_);
3753  for (int j=0; j<curDim_; ++j) {
3754  for (int i=0; i<curDim_; ++i) {
3755  SDMerr(i,j) = subKK(i,j) - SCT::conjugate(subKK(j,i));
3756  }
3757  }
3758  os << " >> Error in KK - KK^H == 0 : " << SDMerr.normFrobenius() << endl;
3759  }
3760 
3761  // Q
3762  if (chk.checkQ) {
3763  for (Array_size_type i=0; i<auxVecs_.size(); ++i) {
3764  MagnitudeType err = orthman_->orthonormError(*auxVecs_[i]);
3765  os << " >> Error in Q[" << i << "]^H M Q[" << i << "] == I : " << err << endl;
3766  for (Array_size_type j=i+1; j<auxVecs_.size(); ++j) {
3767  err = orthman_->orthogError(*auxVecs_[i],*auxVecs_[j]);
3768  os << " >> Error in Q[" << i << "]^H M Q[" << j << "] == 0 : " << err << endl;
3769  }
3770  }
3771  }
3772 
3773  os << endl;
3774 
3775  return os.str();
3776  }
3777 
3778 }} // End of namespace Anasazi
3779 
3780 #endif
3781 
3782 // End of file AnasaziTraceMinBase.hpp
RCP< const MV > getRitzVectors()
Get access to the current Ritz vectors.
int NEV
Number of unconverged eigenvalues.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
RCP< const MV > X
The current eigenvectors.
This class defines the interface required by an eigensolver and status test class to compute solution...
Structure to contain pointers to TraceMinBase state variables.
An implementation of the Anasazi::SortManager that performs a collection of common sorting techniques...
void setBlockSize(int blockSize)
Set the blocksize.
RCP< const MV > V
The current basis.
Declaration of basic traits for the multivector type.
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...
int getBlockSize() const
Get the blocksize used by the iterative solver.
T & get(const std::string &name, T def_value)
An operator of the form [A Y; Y&#39; 0] where A is a sparse matrix and Y a multivector.
RCP< const MV > MopV
The image of V under M, or Teuchos::null if M was not specified.
TraceMinBaseState< ScalarType, MV > getState() const
Get access to the current state of the eigensolver.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
TraceMinBaseOrthoFailure is thrown when the orthogonalization manager is unable to orthogonalize the ...
Virtual base class which defines basic traits for the operator type.
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...
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace. For the trace minimization methods...
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > getRitzRes2Norms()
Get the 2-norms of the residuals.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Basic implementation of the Anasazi::SortManager class.
TraceMinBaseInitFailure is thrown when the TraceMinBase solver is unable to generate an initial itera...
An exception class parent to all Anasazi exceptions.
std::vector< int > getRitzIndex()
Get the index used for extracting individual Ritz vectors from getRitzVectors().
void setStatusTest(RCP< StatusTest< ScalarType, MV, OP > > test)
Set a new StatusTest for the solver.
void initialize()
Initialize the solver with the initial vectors from the eigenproblem or random data.
void currentStatus(std::ostream &os)
This method requests that the solver print out its current status to the given output stream...
bool isOrtho
Whether V has been projected and orthonormalized already.
Pure virtual base class which describes the basic interface to the iterative eigensolver.
Anasazi&#39;s templated, static class providing utilities for the solvers.
bool isInitialized() const
Indicates whether the solver has been initialized or not.
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)
virtual ~TraceMinBase()
Anasazi::TraceMinBase destructor.
Output managers remove the need for the eigensolver to know any information about the required output...
Teuchos::Array< RCP< const MV > > getAuxVecs() const
Get the auxiliary vectors for the solver.
Traits class which defines basic operations on multivectors.
void iterate()
This method performs trace minimization iterations until the status test indicates the need to stop o...
Virtual base class which defines basic traits for the operator type.
int curDim
The current dimension of the solver.
RCP< const MV > KV
The image of V under K.
Anasazi header file which uses auto-configuration information to include necessary C++ headers...
void resize(size_type new_size, const value_type &x=value_type())
RCP< const MV > MX
The image of the current eigenvectors under M, or Teuchos::null if M was not specified.
RCP< StatusTest< ScalarType, MV, OP > > getStatusTest() const
Get the current StatusTest used by the solver.
RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > RV
The current Ritz vectors.
ScalarType largestSafeShift
Largest safe shift.
RCP< const std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > > T
The current Ritz values. This vector is a copy of the internal data.
std::vector< Value< ScalarType > > getRitzValues()
Get the Ritz values for the previous iteration.
RCP< const std::vector< ScalarType > > ritzShifts
Current Ritz shifts.
void sort(std::vector< MagnitudeType > &evals, Teuchos::RCP< std::vector< int > > perm=Teuchos::null, int n=-1) const
Sort real eigenvalues, optionally returning the permutation vector.
RCP< const MV > KX
The image of the current eigenvectors under K.
RCP< const MV > R
The current residual vectors.
RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > KK
The current projected K matrix.
Types and exceptions used within Anasazi solvers and interfaces.
const Eigenproblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the eigenvalue problem.
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 is an abstract base class for the trace minimization eigensolvers.
int getNumIters() const
Get the current iteration count.
Anasazi&#39;s templated pure virtual class for managing the sorting of approximate eigenvalues computed b...
Stores a set of vectors of the form [V; L] where V is a distributed multivector and L is a serialdens...
Common interface of stopping criteria for Anasazi&#39;s solvers.
void setAuxVecs(const Teuchos::Array< RCP< const MV > > &auxvecs)
Set the auxiliary vectors for the solver.
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current eigenvectors and eigenvalues...
The Eigensolver is a templated virtual base class that defines the basic interface that any eigensolv...
void resetNumIters()
Reset the iteration count.
TraceMinBase(const RCP< Eigenproblem< ScalarType, MV, OP > > &problem, const RCP< SortManager< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > > &sorter, const RCP< OutputManager< ScalarType > > &printer, const RCP< StatusTest< ScalarType, MV, OP > > &tester, const RCP< MatOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList &params)
TraceMinBase constructor with eigenproblem, solver utilities, and parameter list of solver options...
int setMatrix(const RCP< SerialDenseMatrix< OrdinalType, ScalarType > > &A)
Class which provides internal utilities for the Anasazi solvers.