Anasazi  Version of the Day
AnasaziGeneralizedDavidson.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Anasazi: Block Eigensolvers Package
5 // Copyright 2004 Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef ANASAZI_GENERALIZED_DAVIDSON_HPP
43 #define ANASAZI_GENERALIZED_DAVIDSON_HPP
44 
51 #include "Teuchos_RCPDecl.hpp"
52 #include "Teuchos_ParameterList.hpp"
53 #include "Teuchos_SerialDenseMatrix.hpp"
54 #include "Teuchos_SerialDenseVector.hpp"
55 #include "Teuchos_Array.hpp"
56 #include "Teuchos_BLAS.hpp"
57 #include "Teuchos_LAPACK.hpp"
58 
59 #include "AnasaziConfigDefs.hpp"
60 #include "AnasaziTypes.hpp"
61 #include "AnasaziEigenproblem.hpp"
62 #include "AnasaziEigensolver.hpp"
63 #include "AnasaziOrthoManager.hpp"
64 #include "AnasaziOutputManager.hpp"
65 #include "AnasaziSortManager.hpp"
66 #include "AnasaziStatusTest.hpp"
67 
68 using Teuchos::RCP;
69 
70 namespace Anasazi {
71 
75 template <class ScalarType, class MV>
78  int curDim;
79 
81  RCP<MV> V;
82 
84  RCP<MV> AV;
85 
87  RCP<MV> BV;
88 
90  RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > VAV;
91 
93  RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > VBV;
94 
96  RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > S;
97 
99  RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > T;
100 
102  RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > Q;
103 
105  RCP< Teuchos::SerialDenseMatrix<int,ScalarType> > Z;
106 
108  std::vector< Value<ScalarType> > eVals;
109 
110  GeneralizedDavidsonState() : curDim(0), V(Teuchos::null), AV(Teuchos::null),
111  BV(Teuchos::null), VAV(Teuchos::null),
112  VBV(Teuchos::null), S(Teuchos::null),
113  T(Teuchos::null), Q(Teuchos::null),
114  Z(Teuchos::null), eVals(0) {}
115 
116 };
117 
118 
139 template <class ScalarType, class MV, class OP>
140 class GeneralizedDavidson : public Eigensolver<ScalarType,MV,OP>
141 {
142  private:
143  // Convenience Typedefs
146  typedef Teuchos::ScalarTraits<ScalarType> ST;
147  typedef typename ST::magnitudeType MagnitudeType;
148  typedef Teuchos::ScalarTraits<MagnitudeType> MT;
149 
150  public:
151 
173  const RCP<SortManager<MagnitudeType> > &sortman,
174  const RCP<OutputManager<ScalarType> > &outputman,
175  const RCP<StatusTest<ScalarType,MV,OP> > &tester,
176  const RCP<OrthoManager<ScalarType,MV> > &orthoman,
177  Teuchos::ParameterList &pl);
178 
182  void iterate();
183 
193  void initialize();
194 
198  void initialize( GeneralizedDavidsonState<ScalarType,MV>& state );
199 
203  int getNumIters() const { return d_iteration; }
204 
208  void resetNumIters() { d_iteration=0; d_opApplies=0; }
209 
213  RCP<const MV> getRitzVectors()
214  {
215  if( !d_ritzVectorsValid )
216  computeRitzVectors();
217  return d_ritzVecs;
218  }
219 
223  std::vector< Value<ScalarType> > getRitzValues();
224 
228  std::vector<int> getRitzIndex()
229  {
230  if( !d_ritzIndexValid )
231  computeRitzIndex();
232  return d_ritzIndex;
233  }
234 
240  std::vector<int> getBlockIndex() const
241  {
242  return d_expansionIndices;
243  }
244 
248  std::vector<MagnitudeType> getResNorms();
249 
253  std::vector<MagnitudeType> getResNorms(int numWanted);
254 
258  std::vector<MagnitudeType> getRes2Norms() { return d_resNorms; }
259 
266  std::vector<MagnitudeType> getRitzRes2Norms() { return d_resNorms; }
267 
271  int getCurSubspaceDim() const { return d_curDim; }
272 
276  int getMaxSubspaceDim() const { return d_maxSubspaceDim; }
277 
281  void setStatusTest( RCP<StatusTest<ScalarType,MV,OP> > tester) { d_tester = tester; }
282 
286  RCP<StatusTest<ScalarType,MV,OP> > getStatusTest() const { return d_tester; }
287 
291  const Eigenproblem<ScalarType,MV,OP> & getProblem() const { return *d_problem; }
292 
296  int getBlockSize() const { return d_expansionSize; }
297 
301  void setBlockSize(int blockSize);
302 
306  void setSize(int blockSize, int maxSubDim);
307 
311  Teuchos::Array< RCP<const MV> > getAuxVecs() const { return d_auxVecs; }
312 
320  void setAuxVecs( const Teuchos::Array< RCP<const MV> > &auxVecs );
321 
325  bool isInitialized() const { return d_initialized; }
326 
330  void currentStatus( std::ostream &myout );
331 
336 
340  void sortProblem( int numWanted );
341 
342  private:
343 
344  // Expand subspace
345  void expandSearchSpace();
346 
347  // Apply Operators
348  void applyOperators();
349 
350  // Update projections
351  void updateProjections();
352 
353  // Solve projected eigenproblem
354  void solveProjectedEigenproblem();
355 
356  // Compute eigenvectors of matrix pair
357  void computeProjectedEigenvectors( Teuchos::SerialDenseMatrix<int,ScalarType> &X );
358 
359  // Scale projected eigenvectors by alpha/beta
360  void scaleEigenvectors( const Teuchos::SerialDenseMatrix<int,ScalarType> &X,
361  Teuchos::SerialDenseMatrix<int,ScalarType> &X_alpha,
362  Teuchos::SerialDenseMatrix<int,ScalarType> &X_beta );
363 
364  // Sort vectors of pairs
365  void sortValues( std::vector<MagnitudeType> &realParts,
366  std::vector<MagnitudeType> &imagParts,
367  std::vector<int> &permVec,
368  int N);
369 
370  // Compute Residual
371  void computeResidual();
372 
373  // Update the current Ritz index vector
374  void computeRitzIndex();
375 
376  // Compute the current Ritz vectors
377  void computeRitzVectors();
378 
379  // Operators
380  RCP<Eigenproblem<ScalarType,MV,OP> > d_problem;
381  Teuchos::ParameterList d_pl;
382  RCP<const OP> d_A;
383  RCP<const OP> d_B;
384  RCP<const OP> d_P;
385  bool d_haveB;
386  bool d_haveP;
387 
388  // Parameters
389  int d_blockSize;
390  int d_maxSubspaceDim;
391  int d_NEV;
392  int d_numToPrint;
393  std::string d_initType;
394  int d_verbosity;
395  bool d_relativeConvergence;
396 
397  // Managers
398  RCP<OutputManager<ScalarType> > d_outputMan;
399  RCP<OrthoManager<ScalarType,MV> > d_orthoMan;
400  RCP<SortManager<MagnitudeType> > d_sortMan;
401  RCP<StatusTest<ScalarType,MV,OP> > d_tester;
402 
403  // Eigenvalues
404  std::vector< Value<ScalarType> > d_eigenvalues;
405 
406  // Residual Vector
407  RCP<MV> d_R;
408  std::vector<MagnitudeType> d_resNorms;
409 
410  // Subspace Vectors
411  RCP<MV> d_V;
412  RCP<MV> d_AV;
413  RCP<MV> d_BV;
414  RCP<MV> d_ritzVecSpace;
415  RCP<MV> d_ritzVecs;
416  Teuchos::Array< RCP<const MV> > d_auxVecs;
417 
418  // Serial Matrices
419  RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > d_VAV;
420  RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > d_VBV;
421  RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > d_S;
422  RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > d_T;
423  RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > d_Q;
424  RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > d_Z;
425 
426  // Arrays for holding Ritz values
427  std::vector<MagnitudeType> d_alphar;
428  std::vector<MagnitudeType> d_alphai;
429  std::vector<MagnitudeType> d_betar;
430  std::vector<int> d_ritzIndex;
431  std::vector<int> d_convergedIndices;
432  std::vector<int> d_expansionIndices;
433 
434  // Current subspace dimension
435  int d_curDim;
436 
437  // How many vectors are to be added to the subspace
438  int d_expansionSize;
439 
440  // Should subspace expansion use leading vectors
441  // (if false, will use leading unconverged vectors)
442  bool d_useLeading;
443 
444  // What should be used for test subspace (V, AV, or BV)
445  std::string d_testSpace;
446 
447  // How many residual vectors are valid
448  int d_residualSize;
449 
450  int d_iteration;
451  int d_opApplies;
452  bool d_initialized;
453  bool d_ritzIndexValid;
454  bool d_ritzVectorsValid;
455 
456 };
457 
458 //---------------------------------------------------------------------------//
459 // Prevent instantiation on complex scalar type
460 //---------------------------------------------------------------------------//
461 template <class MagnitudeType, class MV, class OP>
462 class GeneralizedDavidson<std::complex<MagnitudeType>,MV,OP>
463 {
464  public:
465 
466  typedef std::complex<MagnitudeType> ScalarType;
468  const RCP<Eigenproblem<ScalarType,MV,OP> > &problem,
469  const RCP<SortManager<MagnitudeType> > &sortman,
470  const RCP<OutputManager<ScalarType> > &outputman,
471  const RCP<StatusTest<ScalarType,MV,OP> > &tester,
472  const RCP<OrthoManager<ScalarType,MV> > &orthoman,
473  Teuchos::ParameterList &pl)
474  {
475  // Provide a compile error when attempting to instantiate on complex type
476  MagnitudeType::this_class_is_missing_a_specialization();
477  }
478 };
479 
480 //---------------------------------------------------------------------------//
481 // PUBLIC METHODS
482 //---------------------------------------------------------------------------//
483 
484 //---------------------------------------------------------------------------//
485 // Constructor
486 //---------------------------------------------------------------------------//
487 template <class ScalarType, class MV, class OP>
489  const RCP<Eigenproblem<ScalarType,MV,OP> > &problem,
490  const RCP<SortManager<MagnitudeType> > &sortman,
491  const RCP<OutputManager<ScalarType> > &outputman,
492  const RCP<StatusTest<ScalarType,MV,OP> > &tester,
493  const RCP<OrthoManager<ScalarType,MV> > &orthoman,
494  Teuchos::ParameterList &pl )
495 {
496  TEUCHOS_TEST_FOR_EXCEPTION( problem == Teuchos::null, std::invalid_argument, "No Eigenproblem given to solver." );
497  TEUCHOS_TEST_FOR_EXCEPTION( outputman == Teuchos::null, std::invalid_argument, "No OutputManager given to solver." );
498  TEUCHOS_TEST_FOR_EXCEPTION( orthoman == Teuchos::null, std::invalid_argument, "No OrthoManager given to solver." );
499  TEUCHOS_TEST_FOR_EXCEPTION( sortman == Teuchos::null, std::invalid_argument, "No SortManager given to solver." );
500  TEUCHOS_TEST_FOR_EXCEPTION( tester == Teuchos::null, std::invalid_argument, "No StatusTest given to solver." );
501  TEUCHOS_TEST_FOR_EXCEPTION( !problem->isProblemSet(), std::invalid_argument, "Problem has not been set." );
502 
503  d_problem = problem;
504  d_pl = pl;
505  TEUCHOS_TEST_FOR_EXCEPTION( problem->getA()==Teuchos::null && problem->getOperator()==Teuchos::null,
506  std::invalid_argument, "Either A or Operator must be non-null on Eigenproblem");
507  d_A = problem->getA()!=Teuchos::null ? problem->getA() : problem->getOperator();
508  d_B = problem->getM();
509  d_P = problem->getPrec();
510  d_sortMan = sortman;
511  d_outputMan = outputman;
512  d_tester = tester;
513  d_orthoMan = orthoman;
514 
515  // Pull entries from the ParameterList and Eigenproblem
516  d_NEV = d_problem->getNEV();
517  d_initType = d_pl.get<std::string>("Initial Guess","Random");
518  d_numToPrint = d_pl.get<int>("Print Number of Ritz Values",-1);
519  d_useLeading = d_pl.get<bool>("Use Leading Vectors",false);
520 
521  if( d_B != Teuchos::null )
522  d_haveB = true;
523  else
524  d_haveB = false;
525 
526  if( d_P != Teuchos::null )
527  d_haveP = true;
528  else
529  d_haveP = false;
530 
531  d_testSpace = d_pl.get<std::string>("Test Space","V");
532  TEUCHOS_TEST_FOR_EXCEPTION( d_testSpace!="V" && d_testSpace!="AV" && d_testSpace!="BV", std::invalid_argument,
533  "Anasazi::GeneralizedDavidson: Test Space must be V, AV, or BV" );
534  TEUCHOS_TEST_FOR_EXCEPTION( d_testSpace=="V" ? false : !d_haveB, std::invalid_argument,
535  "Anasazi::GeneralizedDavidson: Test Space must be V for standard eigenvalue problem" );
536 
537  // Allocate space for subspace vectors, projected matrices
538  int blockSize = d_pl.get<int>("Block Size",1);
539  int maxSubDim = d_pl.get<int>("Maximum Subspace Dimension",3*d_NEV*blockSize);
540  d_blockSize = -1;
541  d_maxSubspaceDim = -1;
542  setSize( blockSize, maxSubDim );
543  d_relativeConvergence = d_pl.get<bool>("Relative Convergence Tolerance",false);
544 
545  // Make sure subspace size is consistent with requested eigenvalues
546  TEUCHOS_TEST_FOR_EXCEPTION( d_blockSize <= 0, std::invalid_argument, "Block size must be positive");
547  TEUCHOS_TEST_FOR_EXCEPTION( d_maxSubspaceDim <= 0, std::invalid_argument, "Maximum Subspace Dimension must be positive" );
548  TEUCHOS_TEST_FOR_EXCEPTION( d_problem->getNEV()+2 > pl.get<int>("Maximum Subspace Dimension"),
549  std::invalid_argument, "Maximum Subspace Dimension must be strictly greater than NEV");
550  TEUCHOS_TEST_FOR_EXCEPTION( d_maxSubspaceDim > MVT::GetGlobalLength(*problem->getInitVec()),
551  std::invalid_argument, "Maximum Subspace Dimension cannot exceed problem size");
552 
553 
554  d_curDim = 0;
555  d_iteration = 0;
556  d_opApplies = 0;
557  d_ritzIndexValid = false;
558  d_ritzVectorsValid = false;
559 }
560 
561 
562 //---------------------------------------------------------------------------//
563 // Iterate
564 //---------------------------------------------------------------------------//
565 template <class ScalarType, class MV, class OP>
567 {
568  // Initialize Problem
569  if( !d_initialized )
570  {
571  d_outputMan->stream(Warnings) << "WARNING: GeneralizedDavidson::iterate called without first calling initialize" << std::endl;
572  d_outputMan->stream(Warnings) << " Default initialization will be performed" << std::endl;
573  initialize();
574  }
575 
576  // Print current status
577  if( d_outputMan->isVerbosity(Debug) )
578  {
579  currentStatus( d_outputMan->stream(Debug) );
580  }
581  else if( d_outputMan->isVerbosity(IterationDetails) )
582  {
583  currentStatus( d_outputMan->stream(IterationDetails) );
584  }
585 
586  while( d_tester->getStatus() != Passed && d_curDim+d_expansionSize <= d_maxSubspaceDim )
587  {
588  d_iteration++;
589 
590  expandSearchSpace();
591 
592  applyOperators();
593 
594  updateProjections();
595 
596  solveProjectedEigenproblem();
597 
598  // Make sure the most significant Ritz values are in front
599  // We want the greater of the block size and the number of
600  // requested values, but can't exceed the current dimension
601  int numToSort = std::max(d_blockSize,d_NEV);
602  numToSort = std::min(numToSort,d_curDim);
603  sortProblem( numToSort );
604 
605  computeResidual();
606 
607  // Print current status
608  if( d_outputMan->isVerbosity(Debug) )
609  {
610  currentStatus( d_outputMan->stream(Debug) );
611  }
612  else if( d_outputMan->isVerbosity(IterationDetails) )
613  {
614  currentStatus( d_outputMan->stream(IterationDetails) );
615  }
616  }
617 }
618 
619 //---------------------------------------------------------------------------//
620 // Return the current state struct
621 //---------------------------------------------------------------------------//
622 template <class ScalarType, class MV, class OP>
624 {
626  state.curDim = d_curDim;
627  state.V = d_V;
628  state.AV = d_AV;
629  state.BV = d_BV;
630  state.VAV = d_VAV;
631  state.VBV = d_VBV;
632  state.S = d_S;
633  state.T = d_T;
634  state.Q = d_Q;
635  state.Z = d_Z;
636  state.eVals = getRitzValues();
637  return state;
638 }
639 
640 //---------------------------------------------------------------------------//
641 // Set block size
642 //---------------------------------------------------------------------------//
643 template <class ScalarType, class MV, class OP>
645 {
646  setSize(blockSize,d_maxSubspaceDim);
647 }
648 
649 //---------------------------------------------------------------------------//
650 // Set block size and maximum subspace dimension.
651 //---------------------------------------------------------------------------//
652 template <class ScalarType, class MV, class OP>
653 void GeneralizedDavidson<ScalarType,MV,OP>::setSize(int blockSize, int maxSubDim )
654 {
655  if( blockSize != d_blockSize || maxSubDim != d_maxSubspaceDim )
656  {
657  d_blockSize = blockSize;
658  d_maxSubspaceDim = maxSubDim;
659  d_initialized = false;
660 
661  d_outputMan->stream(Debug) << " >> Anasazi::GeneralizedDavidson: Allocating eigenproblem"
662  << " state with block size of " << d_blockSize
663  << " and maximum subspace dimension of " << d_maxSubspaceDim << std::endl;
664 
665  // Resize arrays for Ritz values
666  d_alphar.resize(d_maxSubspaceDim);
667  d_alphai.resize(d_maxSubspaceDim);
668  d_betar.resize(d_maxSubspaceDim);
669 
670  // Shorten for convenience here
671  int msd = d_maxSubspaceDim;
672 
673  // Temporarily save initialization vector to clone needed vectors
674  RCP<const MV> initVec = d_problem->getInitVec();
675 
676  // Allocate subspace vectors
677  d_V = MVT::Clone(*initVec, msd);
678  d_AV = MVT::Clone(*initVec, msd);
679 
680  // Allocate serial matrices
681  d_VAV = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(msd,msd) );
682  d_S = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(msd,msd) );
683  d_Q = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(msd,msd) );
684  d_Z = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(msd,msd) );
685 
686  // If this is generalized eigenproblem, allocate B components
687  if( d_haveB )
688  {
689  d_BV = MVT::Clone(*initVec, msd);
690  d_VBV = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(msd,msd) );
691  d_T = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(msd,msd) );
692  }
693 
694  /* Allocate space for residual and Ritz vectors
695  * The residual serves two purposes in the Davidson algorithm --
696  * subspace expansion (via the preconditioner) and convergence checking.
697  * We need "Block Size" vectors for subspace expantion and NEV vectors
698  * for convergence checking. Allocate space for max of these, one
699  * extra to avoid splitting conjugate pairs
700  * Allocate one more than "Block Size" to avoid splitting a conjugate pair
701  */
702  d_R = MVT::Clone(*initVec,std::max(d_blockSize,d_NEV)+1);
703  d_ritzVecSpace = MVT::Clone(*initVec,std::max(d_blockSize,d_NEV)+1);
704  }
705 }
706 
707 //---------------------------------------------------------------------------//
708 /*
709  * Initialize the eigenvalue problem
710  *
711  * Anything on the state that is not null is assumed to be valid.
712  * Anything not present on the state will be generated
713  * Very limited error checking can be performed to ensure the validity of
714  * state components (e.g. we cannot verify that state.AV actually corresponds
715  * to A*state.V), so this function should be used carefully.
716  */
717 //---------------------------------------------------------------------------//
718 template <class ScalarType, class MV, class OP>
720 {
721  // If state has nonzero dimension, we initialize from that, otherwise
722  // we'll pick d_blockSize vectors to start with
723  d_curDim = (state.curDim > 0 ? state.curDim : d_blockSize );
724 
725  d_outputMan->stream(Debug) << " >> Anasazi::GeneralizedDavidson: Initializing state"
726  << " with subspace dimension " << d_curDim << std::endl;
727 
728  // Index for 1st block_size vectors
729  std::vector<int> initInds(d_curDim);
730  for( int i=0; i<d_curDim; ++i )
731  initInds[i] = i;
732 
733  // View of vectors that need to be initialized
734  RCP<MV> V1 = MVT::CloneViewNonConst(*d_V,initInds);
735 
736  // If state's dimension is large enough, use state.V to initialize
737  bool reset_V = false;;
738  if( state.curDim > 0 && state.V != Teuchos::null && MVT::GetNumberVecs(*state.V) >= d_curDim )
739  {
740  if( state.V != d_V )
741  MVT::SetBlock(*state.V,initInds,*V1);
742  }
743  // If there aren't enough vectors in problem->getInitVec() or the user specifically
744  // wants to use random data, set V to random
745  else if( MVT::GetNumberVecs(*d_problem->getInitVec()) < d_blockSize || d_initType == "Random" )
746  {
747  MVT::MvRandom(*V1);
748  reset_V = true;
749  }
750  // Use vectors in problem->getInitVec()
751  else
752  {
753  RCP<const MV> initVec = MVT::CloneView(*d_problem->getInitVec(),initInds);
754  MVT::SetBlock(*initVec,initInds,*V1);
755  reset_V = true;
756  }
757 
758  // If we reset V, it needs to be orthonormalized
759  if( reset_V )
760  {
761  int rank = d_orthoMan->projectAndNormalize( *V1, d_auxVecs );
762  TEUCHOS_TEST_FOR_EXCEPTION( rank < d_blockSize, std::logic_error,
763  "Anasazi::GeneralizedDavidson::initialize(): Error generating initial orthonormal basis" );
764  }
765 
766  if( d_outputMan->isVerbosity(Debug) )
767  {
768  d_outputMan->stream(Debug) << " >> Anasazi::GeneralizedDavidson: Error in V^T V == I: "
769  << d_orthoMan->orthonormError( *V1 ) << std::endl;
770  }
771 
772  // Now process AV
773  RCP<MV> AV1 = MVT::CloneViewNonConst(*d_AV,initInds);
774 
775  // If AV in the state is valid and of appropriate size, use it
776  // We have no way to check that AV is actually A*V
777  if( !reset_V && state.AV != Teuchos::null && MVT::GetNumberVecs(*state.AV) >= d_curDim )
778  {
779  if( state.AV != d_AV )
780  MVT::SetBlock(*state.AV,initInds,*AV1);
781  }
782  // Otherwise apply A to V
783  else
784  {
785  OPT::Apply( *d_A, *V1, *AV1 );
786  d_opApplies += MVT::GetNumberVecs( *V1 );
787  }
788 
789  // Views of matrix to be updated
790  Teuchos::SerialDenseMatrix<int,ScalarType> VAV1( Teuchos::View, *d_VAV, d_curDim, d_curDim );
791 
792  // If the state has a valid VAV, use it
793  if( !reset_V && state.VAV != Teuchos::null && state.VAV->numRows() >= d_curDim && state.VAV->numCols() >= d_curDim )
794  {
795  if( state.VAV != d_VAV )
796  {
797  Teuchos::SerialDenseMatrix<int,ScalarType> state_VAV( Teuchos::View, *state.VAV, d_curDim, d_curDim );
798  VAV1.assign( state_VAV );
799  }
800  }
801  // Otherwise compute VAV from V,AV
802  else
803  {
804  if( d_testSpace == "V" )
805  {
806  MVT::MvTransMv( ST::one(), *V1, *AV1, VAV1 );
807  }
808  else if( d_testSpace == "AV" )
809  {
810  MVT::MvTransMv( ST::one(), *AV1, *AV1, VAV1 );
811  }
812  else if( d_testSpace == "BV" )
813  {
814  RCP<MV> BV1 = MVT::CloneViewNonConst(*d_BV,initInds);
815  MVT::MvTransMv( ST::one(), *BV1, *AV1, VAV1 );
816  }
817  }
818 
819  // Process BV if we have it
820  if( d_haveB )
821  {
822  RCP<MV> BV1 = MVT::CloneViewNonConst(*d_BV,initInds);
823 
824  // If BV in the state is valid and of appropriate size, use it
825  // We have no way to check that BV is actually B*V
826  if( !reset_V && state.BV != Teuchos::null && MVT::GetNumberVecs(*state.BV) >= d_curDim )
827  {
828  if( state.BV != d_BV )
829  MVT::SetBlock(*state.BV,initInds,*BV1);
830  }
831  // Otherwise apply B to V
832  else
833  {
834  OPT::Apply( *d_B, *V1, *BV1 );
835  }
836 
837  // Views of matrix to be updated
838  Teuchos::SerialDenseMatrix<int,ScalarType> VBV1( Teuchos::View, *d_VBV, d_curDim, d_curDim );
839 
840  // If the state has a valid VBV, use it
841  if( !reset_V && state.VBV != Teuchos::null && state.VBV->numRows() >= d_curDim && state.VBV->numCols() >= d_curDim )
842  {
843  if( state.VBV != d_VBV )
844  {
845  Teuchos::SerialDenseMatrix<int,ScalarType> state_VBV( Teuchos::View, *state.VBV, d_curDim, d_curDim );
846  VBV1.assign( state_VBV );
847  }
848  }
849  // Otherwise compute VBV from V,BV
850  else
851  {
852  if( d_testSpace == "V" )
853  {
854  MVT::MvTransMv( ST::one(), *V1, *BV1, VBV1 );
855  }
856  else if( d_testSpace == "AV" )
857  {
858  MVT::MvTransMv( ST::one(), *AV1, *BV1, VBV1 );
859  }
860  else if( d_testSpace == "BV" )
861  {
862  MVT::MvTransMv( ST::one(), *BV1, *BV1, VBV1 );
863  }
864  }
865  }
866 
867  // Update Ritz values
868  solveProjectedEigenproblem();
869 
870  // Sort
871  int numToSort = std::max(d_blockSize,d_NEV);
872  numToSort = std::min(numToSort,d_curDim);
873  sortProblem( numToSort );
874 
875  // Get valid residual
876  computeResidual();
877 
878  // Set solver to initialized
879  d_initialized = true;
880 }
881 
882 //---------------------------------------------------------------------------//
883 // Initialize the eigenvalue problem with empty state
884 //---------------------------------------------------------------------------//
885 template <class ScalarType, class MV, class OP>
887 {
889  initialize( empty );
890 }
891 
892 //---------------------------------------------------------------------------//
893 // Get current residual norms
894 //---------------------------------------------------------------------------//
895 template <class ScalarType, class MV, class OP>
896 std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
898 {
899  return getResNorms(d_residualSize);
900 }
901 
902 //---------------------------------------------------------------------------//
903 // Get current residual norms
904 //---------------------------------------------------------------------------//
905 template <class ScalarType, class MV, class OP>
906 std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
908 {
909  std::vector<int> resIndices(numWanted);
910  for( int i=0; i<numWanted; ++i )
911  resIndices[i]=i;
912 
913  RCP<const MV> R_checked = MVT::CloneView( *d_R, resIndices );
914 
915  std::vector<MagnitudeType> resNorms;
916  d_orthoMan->norm( *R_checked, resNorms );
917 
918  return resNorms;
919 }
920 
921 //---------------------------------------------------------------------------//
922 // Get current Ritz values
923 //---------------------------------------------------------------------------//
924 template <class ScalarType, class MV, class OP>
926 {
927  std::vector< Value<ScalarType> > ritzValues;
928  for( int ival=0; ival<d_curDim; ++ival )
929  {
930  Value<ScalarType> thisVal;
931  thisVal.realpart = d_alphar[ival] / d_betar[ival];
932  if( d_betar[ival] != MT::zero() )
933  thisVal.imagpart = d_alphai[ival] / d_betar[ival];
934  else
935  thisVal.imagpart = MT::zero();
936 
937  ritzValues.push_back( thisVal );
938  }
939 
940  return ritzValues;
941 }
942 
943 //---------------------------------------------------------------------------//
944 /*
945  * Set auxilliary vectors
946  *
947  * Manually setting the auxilliary vectors invalidates the current state
948  * of the solver. Reuse of any components of the solver requires extracting
949  * the state, orthogonalizing V against the aux vecs and reinitializing.
950  */
951 //---------------------------------------------------------------------------//
952 template <class ScalarType, class MV, class OP>
954  const Teuchos::Array< RCP<const MV> > &auxVecs )
955 {
956  d_auxVecs = auxVecs;
957 
958  // Set state to uninitialized if any vectors were set here
959  typename Teuchos::Array< RCP<const MV> >::const_iterator arrItr;
960  int numAuxVecs=0;
961  for( arrItr=auxVecs.begin(); arrItr!=auxVecs.end(); ++arrItr )
962  {
963  numAuxVecs += MVT::GetNumberVecs( *(*arrItr) );
964  }
965  if( numAuxVecs > 0 )
966  d_initialized = false;
967 }
968 
969 //---------------------------------------------------------------------------//
970 // Reorder Schur form, bringing wanted values to front
971 //---------------------------------------------------------------------------//
972 template <class ScalarType, class MV, class OP>
974 {
975  // Get permutation vector
976  std::vector<MagnitudeType> realRitz(d_curDim), imagRitz(d_curDim);
977  std::vector< Value<ScalarType> > ritzVals = getRitzValues();
978  for( int i=0; i<d_curDim; ++i )
979  {
980  realRitz[i] = ritzVals[i].realpart;
981  imagRitz[i] = ritzVals[i].imagpart;
982  }
983 
984  std::vector<int> permVec;
985  sortValues( realRitz, imagRitz, permVec, d_curDim );
986 
987  std::vector<int> sel(d_curDim,0);
988  for( int ii=0; ii<numWanted; ++ii )
989  sel[ permVec[ii] ]=1;
990 
991  if( d_haveB )
992  {
993  int ijob = 0; // reorder only, no condition number estimates
994  int wantq = 1; // keep left Schur vectors
995  int wantz = 1; // keep right Schur vectors
996  int work_size=10*d_maxSubspaceDim+16;
997  std::vector<ScalarType> work(work_size);
998  int sdim = 0;
999  int iwork_size = 1;
1000  int iwork;
1001  int info = 0;
1002 
1003  Teuchos::LAPACK<int,ScalarType> lapack;
1004  lapack.TGSEN( ijob, wantq, wantz, &sel[0], d_curDim, d_S->values(), d_S->stride(), d_T->values(), d_T->stride(),
1005  &d_alphar[0], &d_alphai[0], &d_betar[0], d_Q->values(), d_Q->stride(), d_Z->values(), d_Z->stride(),
1006  &sdim, 0, 0, 0, &work[0], work_size, &iwork, iwork_size, &info );
1007 
1008  d_ritzIndexValid = false;
1009  d_ritzVectorsValid = false;
1010 
1011  std::stringstream ss;
1012  ss << "Anasazi::GeneralizedDavidson: TGSEN returned error code " << info << std::endl;
1013  TEUCHOS_TEST_FOR_EXCEPTION( info<0, std::runtime_error, ss.str() );
1014  if( info > 0 )
1015  {
1016  // Only issue a warning for positive error code, this usually indicates
1017  // that the system has not been fully reordered, presumably due to ill-conditioning.
1018  // This is usually not detrimental to the calculation.
1019  d_outputMan->stream(Warnings) << "WARNING: " << ss.str() << std::endl;
1020  d_outputMan->stream(Warnings) << " Problem may not be correctly sorted" << std::endl;
1021  }
1022  }
1023  else {
1024  char getCondNum = 'N'; // no condition number estimates
1025  char getQ = 'V'; // keep Schur vectors
1026  int subDim = 0;
1027  int work_size = d_curDim;
1028  std::vector<ScalarType> work(work_size);
1029  int iwork_size = 1;
1030  int iwork = 0;
1031  int info = 0;
1032 
1033  Teuchos::LAPACK<int,ScalarType> lapack;
1034  lapack.TRSEN (getCondNum, getQ, &sel[0], d_curDim, d_S->values (),
1035  d_S->stride (), d_Z->values (), d_Z->stride (),
1036  &d_alphar[0], &d_alphai[0], &subDim, 0, 0, &work[0],
1037  work_size, &iwork, iwork_size, &info );
1038 
1039  std::fill( d_betar.begin(), d_betar.end(), ST::one() );
1040 
1041  d_ritzIndexValid = false;
1042  d_ritzVectorsValid = false;
1043 
1044  TEUCHOS_TEST_FOR_EXCEPTION(
1045  info < 0, std::runtime_error, "Anasazi::GeneralizedDavidson::"
1046  "sortProblem: LAPACK routine TRSEN returned error code INFO = "
1047  << info << " < 0. This indicates that argument " << -info
1048  << " (counting starts with one) of TRSEN had an illegal value.");
1049 
1050  // LAPACK's documentation suggests that this should only happen
1051  // in the real-arithmetic case, because I only see INFO == 1
1052  // possible for DTRSEN, not for ZTRSEN. Nevertheless, it's
1053  // harmless to check regardless.
1054  TEUCHOS_TEST_FOR_EXCEPTION(
1055  info == 1, std::runtime_error, "Anasazi::GeneralizedDavidson::"
1056  "sortProblem: LAPACK routine TRSEN returned error code INFO = 1. "
1057  "This indicates that the reordering failed because some eigenvalues "
1058  "are too close to separate (the problem is very ill-conditioned).");
1059 
1060  TEUCHOS_TEST_FOR_EXCEPTION(
1061  info > 1, std::logic_error, "Anasazi::GeneralizedDavidson::"
1062  "sortProblem: LAPACK routine TRSEN returned error code INFO = "
1063  << info << " > 1. This should not be possible. It may indicate an "
1064  "error either in LAPACK itself, or in Teuchos' LAPACK wrapper.");
1065  }
1066 }
1067 
1068 
1069 //---------------------------------------------------------------------------//
1070 // PRIVATE IMPLEMENTATION
1071 //---------------------------------------------------------------------------//
1072 
1073 //---------------------------------------------------------------------------//
1074 // Expand subspace using preconditioner and orthogonalize
1075 //---------------------------------------------------------------------------//
1076 template <class ScalarType, class MV, class OP>
1078 {
1079  // Get indices into relevant portion of residual and
1080  // location to be added to search space
1081  std::vector<int> newIndices(d_expansionSize);
1082  for( int i=0; i<d_expansionSize; ++i )
1083  {
1084  newIndices[i] = d_curDim+i;
1085  }
1086 
1087  // Get indices into pre-existing search space
1088  std::vector<int> curIndices(d_curDim);
1089  for( int i=0; i<d_curDim; ++i )
1090  curIndices[i] = i;
1091 
1092  // Get View of vectors
1093  RCP<MV> V_new = MVT::CloneViewNonConst( *d_V, newIndices);
1094  RCP<const MV> V_cur = MVT::CloneView( *d_V, curIndices);
1095  RCP<const MV> R_active = MVT::CloneView( *d_R, d_expansionIndices);
1096 
1097  if( d_haveP )
1098  {
1099  // Apply Preconditioner to Residual
1100  OPT::Apply( *d_P, *R_active, *V_new );
1101  }
1102  else
1103  {
1104  // Just copy the residual
1105  MVT::SetBlock( *R_active, newIndices, *d_V );
1106  }
1107 
1108  // Normalize new vector against existing vectors in V plus auxVecs
1109  Teuchos::Array< RCP<const MV> > against = d_auxVecs;
1110  against.push_back( V_cur );
1111  int rank = d_orthoMan->projectAndNormalize(*V_new,against);
1112 
1113  if( d_outputMan->isVerbosity(Debug) )
1114  {
1115  std::vector<int> allIndices(d_curDim+d_expansionSize);
1116  for( int i=0; i<d_curDim+d_expansionSize; ++i )
1117  allIndices[i]=i;
1118 
1119  RCP<const MV> V_all = MVT::CloneView( *d_V, allIndices );
1120 
1121  d_outputMan->stream(Debug) << " >> Anasazi::GeneralizedDavidson: Error in V^T V == I: "
1122  << d_orthoMan->orthonormError( *V_all ) << std::endl;
1123  }
1124 
1125  TEUCHOS_TEST_FOR_EXCEPTION( rank != d_expansionSize, std::runtime_error,
1126  "Anasazi::GeneralizedDavidson::ExpandSearchSpace(): Orthonormalization of new vectors failed" );
1127 
1128 }
1129 
1130 //---------------------------------------------------------------------------//
1131 // Apply operators
1132 //---------------------------------------------------------------------------//
1133 template <class ScalarType, class MV, class OP>
1135 {
1136  // Get indices for different components
1137  std::vector<int> newIndices(d_expansionSize);
1138  for( int i=0; i<d_expansionSize; ++i )
1139  newIndices[i] = d_curDim+i;
1140 
1141  // Get Views
1142  RCP<const MV> V_new = MVT::CloneView( *d_V, newIndices );
1143  RCP<MV> AV_new = MVT::CloneViewNonConst( *d_AV, newIndices );
1144 
1145  // Multiply by A
1146  OPT::Apply( *d_A, *V_new, *AV_new );
1147  d_opApplies += MVT::GetNumberVecs( *V_new );
1148 
1149  // Multiply by B
1150  if( d_haveB )
1151  {
1152  RCP<MV> BV_new = MVT::CloneViewNonConst( *d_BV, newIndices );
1153  OPT::Apply( *d_B, *V_new, *BV_new );
1154  }
1155 }
1156 
1157 //---------------------------------------------------------------------------//
1158 // Update projected matrices.
1159 //---------------------------------------------------------------------------//
1160 template <class ScalarType, class MV, class OP>
1162 {
1163  // Get indices for different components
1164  std::vector<int> newIndices(d_expansionSize);
1165  for( int i=0; i<d_expansionSize; ++i )
1166  newIndices[i] = d_curDim+i;
1167 
1168  std::vector<int> curIndices(d_curDim);
1169  for( int i=0; i<d_curDim; ++i )
1170  curIndices[i] = i;
1171 
1172  std::vector<int> allIndices(d_curDim+d_expansionSize);
1173  for( int i=0; i<d_curDim+d_expansionSize; ++i )
1174  allIndices[i] = i;
1175 
1176  // Test subspace can be V, AV, or BV
1177  RCP<const MV> W_new, W_all;
1178  if( d_testSpace == "V" )
1179  {
1180  W_new = MVT::CloneView(*d_V, newIndices );
1181  W_all = MVT::CloneView(*d_V, allIndices );
1182  }
1183  else if( d_testSpace == "AV" )
1184  {
1185  W_new = MVT::CloneView(*d_AV, newIndices );
1186  W_all = MVT::CloneView(*d_AV, allIndices );
1187  }
1188  else if( d_testSpace == "BV" )
1189  {
1190  W_new = MVT::CloneView(*d_BV, newIndices );
1191  W_all = MVT::CloneView(*d_BV, allIndices );
1192  }
1193 
1194  // Get views of AV
1195  RCP<const MV> AV_new = MVT::CloneView(*d_AV, newIndices);
1196  RCP<const MV> AV_current = MVT::CloneView(*d_AV, curIndices);
1197 
1198  // Last block_size rows of VAV (minus final entry)
1199  Teuchos::SerialDenseMatrix<int,ScalarType> VAV_lastrow( Teuchos::View, *d_VAV, d_expansionSize, d_curDim, d_curDim, 0 );
1200  MVT::MvTransMv( ST::one(), *W_new, *AV_current, VAV_lastrow );
1201 
1202  // Last block_size columns of VAV
1203  Teuchos::SerialDenseMatrix<int,ScalarType> VAV_lastcol( Teuchos::View, *d_VAV, d_curDim+d_expansionSize, d_expansionSize, 0, d_curDim );
1204  MVT::MvTransMv( ST::one(), *W_all, *AV_new, VAV_lastcol );
1205 
1206  if( d_haveB )
1207  {
1208  // Get views of BV
1209  RCP<const MV> BV_new = MVT::CloneView(*d_BV, newIndices);
1210  RCP<const MV> BV_current = MVT::CloneView(*d_BV, curIndices);
1211 
1212  // Last block_size rows of VBV (minus final entry)
1213  Teuchos::SerialDenseMatrix<int,ScalarType> VBV_lastrow( Teuchos::View, *d_VBV, d_expansionSize, d_curDim, d_curDim, 0 );
1214  MVT::MvTransMv( ST::one(), *W_new, *BV_current, VBV_lastrow );
1215 
1216  // Last block_size columns of VBV
1217  Teuchos::SerialDenseMatrix<int,ScalarType> VBV_lastcol( Teuchos::View, *d_VBV, d_curDim+d_expansionSize, d_expansionSize, 0, d_curDim );
1218  MVT::MvTransMv( ST::one(), *W_all, *BV_new, VBV_lastcol );
1219  }
1220 
1221  // All bases are expanded, increase current subspace dimension
1222  d_curDim += d_expansionSize;
1223 
1224  d_ritzIndexValid = false;
1225  d_ritzVectorsValid = false;
1226 }
1227 
1228 //---------------------------------------------------------------------------//
1229 // Solve low dimensional eigenproblem using LAPACK
1230 //---------------------------------------------------------------------------//
1231 template <class ScalarType, class MV, class OP>
1233 {
1234  if( d_haveB )
1235  {
1236  // VAV and VBV need to stay unchanged, GGES will overwrite
1237  // S and T with the triangular matrices from the generalized
1238  // Schur form
1239  d_S->assign(*d_VAV);
1240  d_T->assign(*d_VBV);
1241 
1242  // Get QZ Decomposition of Projected Problem
1243  char leftVecs = 'V'; // compute left vectors
1244  char rightVecs = 'V'; // compute right vectors
1245  char sortVals = 'N'; // don't sort
1246  int sdim;
1247  // int work_size = 10*d_curDim+16;
1248  Teuchos::LAPACK<int,ScalarType> lapack;
1249  int info;
1250  // workspace query
1251  int work_size = -1;
1252  std::vector<ScalarType> work(1);
1253  lapack.GGES( leftVecs, rightVecs, sortVals, NULL, d_curDim, d_S->values(), d_S->stride(),
1254  d_T->values(), d_T->stride(), &sdim, &d_alphar[0], &d_alphai[0], &d_betar[0],
1255  d_Q->values(), d_Q->stride(), d_Z->values(), d_Z->stride(), &work[0], work_size, 0, &info );
1256  // actual call
1257  work_size = work[0];
1258  work.resize(work_size);
1259  lapack.GGES( leftVecs, rightVecs, sortVals, NULL, d_curDim, d_S->values(), d_S->stride(),
1260  d_T->values(), d_T->stride(), &sdim, &d_alphar[0], &d_alphai[0], &d_betar[0],
1261  d_Q->values(), d_Q->stride(), d_Z->values(), d_Z->stride(), &work[0], work_size, 0, &info );
1262 
1263  d_ritzIndexValid = false;
1264  d_ritzVectorsValid = false;
1265 
1266  std::stringstream ss;
1267  ss << "Anasazi::GeneralizedDavidson: GGES returned error code " << info << std::endl;
1268  TEUCHOS_TEST_FOR_EXCEPTION( info!=0, std::runtime_error, ss.str() );
1269  }
1270  else
1271  {
1272  // VAV needs to stay unchanged, GGES will overwrite
1273  // S with the triangular matrix from the Schur form
1274  d_S->assign(*d_VAV);
1275 
1276  // Get QR Decomposition of Projected Problem
1277  char vecs = 'V'; // compute Schur vectors
1278  int sdim;
1279  int work_size = 3*d_curDim;
1280  std::vector<ScalarType> work(work_size);
1281  int info;
1282 
1283  Teuchos::LAPACK<int,ScalarType> lapack;
1284  lapack.GEES( vecs, d_curDim, d_S->values(), d_S->stride(), &sdim, &d_alphar[0], &d_alphai[0],
1285  d_Z->values(), d_Z->stride(), &work[0], work_size, 0, 0, &info);
1286 
1287  std::fill( d_betar.begin(), d_betar.end(), ST::one() );
1288 
1289  d_ritzIndexValid = false;
1290  d_ritzVectorsValid = false;
1291 
1292  std::stringstream ss;
1293  ss << "Anasazi::GeneralizedDavidson: GEES returned error code " << info << std::endl;
1294  TEUCHOS_TEST_FOR_EXCEPTION( info!=0, std::runtime_error, ss.str() );
1295  }
1296 }
1297 
1298 //---------------------------------------------------------------------------//
1299 /*
1300  * Get index vector into current Ritz values/vectors
1301  *
1302  * The current ordering of d_alphar, d_alphai, d_betar will be used.
1303  * Reordering those vectors will invalidate the vector returned here.
1304  */
1305 //---------------------------------------------------------------------------//
1306 template <class ScalarType, class MV, class OP>
1308 {
1309  if( d_ritzIndexValid )
1310  return;
1311 
1312  d_ritzIndex.resize( d_curDim );
1313  int i=0;
1314  while( i < d_curDim )
1315  {
1316  if( d_alphai[i] == ST::zero() )
1317  {
1318  d_ritzIndex[i] = 0;
1319  i++;
1320  }
1321  else
1322  {
1323  d_ritzIndex[i] = 1;
1324  d_ritzIndex[i+1] = -1;
1325  i+=2;
1326  }
1327  }
1328  d_ritzIndexValid = true;
1329 }
1330 
1331 //---------------------------------------------------------------------------//
1332 /*
1333  * Compute current Ritz vectors
1334  *
1335  * The current ordering of d_alphar, d_alphai, d_betar will be used.
1336  * Reordering those vectors will invalidate the vector returned here.
1337  */
1338 //---------------------------------------------------------------------------//
1339 template <class ScalarType, class MV, class OP>
1341 {
1342  if( d_ritzVectorsValid )
1343  return;
1344 
1345  // Make Ritz indices current
1346  computeRitzIndex();
1347 
1348  // Get indices of converged vector
1349  std::vector<int> checkedIndices(d_residualSize);
1350  for( int ii=0; ii<d_residualSize; ++ii )
1351  checkedIndices[ii] = ii;
1352 
1353  // Get eigenvectors of projected system
1354  Teuchos::SerialDenseMatrix<int,ScalarType> X(Teuchos::Copy,*d_Z,d_curDim,d_curDim);
1355  computeProjectedEigenvectors( X );
1356 
1357  // Get view of wanted vectors
1358  Teuchos::SerialDenseMatrix<int,ScalarType> X_wanted(Teuchos::View,X,d_curDim,d_residualSize);
1359 
1360  // Get views of relevant portion of V, evecs
1361  d_ritzVecs = MVT::CloneViewNonConst( *d_ritzVecSpace, checkedIndices );
1362 
1363  std::vector<int> curIndices(d_curDim);
1364  for( int i=0; i<d_curDim; ++i )
1365  curIndices[i] = i;
1366 
1367  RCP<const MV> V_current = MVT::CloneView( *d_V, curIndices );
1368 
1369  // Now form Ritz vector
1370  MVT::MvTimesMatAddMv(ST::one(),*V_current,X_wanted,ST::zero(),*d_ritzVecs);
1371 
1372  // Normalize vectors, conjugate pairs get normalized together
1373  std::vector<MagnitudeType> scale(d_residualSize);
1374  MVT::MvNorm( *d_ritzVecs, scale );
1375  Teuchos::LAPACK<int,ScalarType> lapack;
1376  for( int i=0; i<d_residualSize; ++i )
1377  {
1378  if( d_ritzIndex[i] == 0 )
1379  {
1380  scale[i] = 1.0/scale[i];
1381  }
1382  else if( d_ritzIndex[i] == 1 )
1383  {
1384  MagnitudeType nrm = lapack.LAPY2(scale[i],scale[i+1]);
1385  scale[i] = 1.0/nrm;
1386  scale[i+1] = 1.0/nrm;
1387  }
1388  }
1389  MVT::MvScale( *d_ritzVecs, scale );
1390 
1391  d_ritzVectorsValid = true;
1392 
1393 }
1394 
1395 //---------------------------------------------------------------------------//
1396 // Use sort manager to sort generalized eigenvalues
1397 //---------------------------------------------------------------------------//
1398 template <class ScalarType, class MV, class OP>
1399 void GeneralizedDavidson<ScalarType,MV,OP>::sortValues( std::vector<MagnitudeType> &realParts,
1400  std::vector<MagnitudeType> &imagParts,
1401  std::vector<int> &permVec,
1402  int N)
1403 {
1404  permVec.resize(N);
1405 
1406  TEUCHOS_TEST_FOR_EXCEPTION( (int) realParts.size()<N, std::runtime_error,
1407  "Anasazi::GeneralizedDavidson::SortValues: Number of requested sorted values greater than vector length." );
1408  TEUCHOS_TEST_FOR_EXCEPTION( (int) imagParts.size()<N, std::runtime_error,
1409  "Anasazi::GeneralizedDavidson::SortValues: Number of requested sorted values greater than vector length." );
1410 
1411  RCP< std::vector<int> > rcpPermVec = Teuchos::rcpFromRef(permVec);
1412 
1413  d_sortMan->sort( realParts, imagParts, rcpPermVec, N );
1414 
1415  d_ritzIndexValid = false;
1416  d_ritzVectorsValid = false;
1417 }
1418 
1419 //---------------------------------------------------------------------------//
1420 /*
1421  * Compute (right) scaled eigenvectors of a pair of dense matrices
1422  *
1423  * This routine computes the eigenvectors for the generalized eigenvalue
1424  * problem \f$ \beta A x = \alpha B x \f$. The input matrices are the upper
1425  * quasi-triangular matrices S and T from a real QZ decomposition,
1426  * the routine dtgevc will back-calculate the eigenvectors of the original
1427  * pencil (A,B) using the orthogonal matrices Q and Z.
1428  */
1429 //---------------------------------------------------------------------------//
1430 template <class ScalarType, class MV, class OP>
1432  Teuchos::SerialDenseMatrix<int,ScalarType> &X )
1433 {
1434  int N = X.numRows();
1435  if( d_haveB )
1436  {
1437  Teuchos::SerialDenseMatrix<int,ScalarType> S(Teuchos::Copy, *d_S, N, N);
1438  Teuchos::SerialDenseMatrix<int,ScalarType> T(Teuchos::Copy, *d_T, N, N);
1439  Teuchos::SerialDenseMatrix<int,ScalarType> VL(Teuchos::Copy, *d_Q, N, N);
1440 
1441  char whichVecs = 'R'; // only need right eigenvectors
1442  char howMany = 'B'; // back-compute eigenvectors of original A,B (we have Schur decomposition here)
1443  int work_size = 6*d_maxSubspaceDim;
1444  std::vector<ScalarType> work(work_size,ST::zero());
1445  int info;
1446  int M;
1447 
1448  Teuchos::LAPACK<int,ScalarType> lapack;
1449  lapack.TGEVC( whichVecs, howMany, 0, N, S.values(), S.stride(), T.values(), T.stride(),
1450  VL.values(), VL.stride(), X.values(), X.stride(), N, &M, &work[0], &info );
1451 
1452  std::stringstream ss;
1453  ss << "Anasazi::GeneralizedDavidson: TGEVC returned error code " << info << std::endl;
1454  TEUCHOS_TEST_FOR_EXCEPTION( info!=0, std::runtime_error, ss.str() );
1455  }
1456  else
1457  {
1458  Teuchos::SerialDenseMatrix<int,ScalarType> S(Teuchos::Copy, *d_S, N, N);
1459  Teuchos::SerialDenseMatrix<int,ScalarType> VL(Teuchos::Copy, *d_Z, N, N);
1460 
1461  char whichVecs = 'R'; // only need right eigenvectors
1462  char howMany = 'B'; // back-compute eigenvectors of original A (we have Schur decomposition here)
1463  int sel = 0;
1464  std::vector<ScalarType> work(3*N);
1465  int m;
1466  int info;
1467 
1468  Teuchos::LAPACK<int,ScalarType> lapack;
1469 
1470  lapack.TREVC( whichVecs, howMany, &sel, N, S.values(), S.stride(), VL.values(), VL.stride(),
1471  X.values(), X.stride(), N, &m, &work[0], &info );
1472 
1473  std::stringstream ss;
1474  ss << "Anasazi::GeneralizedDavidson: TREVC returned error code " << info << std::endl;
1475  TEUCHOS_TEST_FOR_EXCEPTION( info!=0, std::runtime_error, ss.str() );
1476  }
1477 }
1478 
1479 //---------------------------------------------------------------------------//
1480 // Scale eigenvectors by quasi-diagonal matrices alpha and beta
1481 //---------------------------------------------------------------------------//
1482 template <class ScalarType, class MV, class OP>
1484  const Teuchos::SerialDenseMatrix<int,ScalarType> &X,
1485  Teuchos::SerialDenseMatrix<int,ScalarType> &X_alpha,
1486  Teuchos::SerialDenseMatrix<int,ScalarType> &X_beta )
1487 {
1488  int Nr = X.numRows();
1489  int Nc = X.numCols();
1490 
1491  TEUCHOS_TEST_FOR_EXCEPTION( Nr>d_curDim, std::logic_error,
1492  "Anasazi::GeneralizedDavidson::ScaleEigenvectors: Matrix size exceeds current dimension");
1493  TEUCHOS_TEST_FOR_EXCEPTION( Nc>d_curDim, std::logic_error,
1494  "Anasazi::GeneralizedDavidson::ScaleEigenvectors: Matrix size exceeds current dimension");
1495  TEUCHOS_TEST_FOR_EXCEPTION( X_alpha.numRows()!=Nr, std::logic_error,
1496  "Anasazi::GeneralizedDavidson::ScaleEigenvectors: number of rows in Xalpha does not match X");
1497  TEUCHOS_TEST_FOR_EXCEPTION( X_alpha.numCols()!=Nc, std::logic_error,
1498  "Anasazi::GeneralizedDavidson::ScaleEigenvectors: number of cols in Xalpha does not match X");
1499  TEUCHOS_TEST_FOR_EXCEPTION( X_beta.numRows()!=Nr, std::logic_error,
1500  "Anasazi::GeneralizedDavidson::ScaleEigenvectors: number of rows in Xbeta does not match X");
1501  TEUCHOS_TEST_FOR_EXCEPTION( X_beta.numCols()!=Nc, std::logic_error,
1502  "Anasazi::GeneralizedDavidson::ScaleEigenvectors: number of cols in Xbeta does not match X");
1503 
1504  // Now form quasi-diagonal matrices
1505  // containing alpha and beta
1506  Teuchos::SerialDenseMatrix<int,ScalarType> Alpha(Nc,Nc,true);
1507  Teuchos::SerialDenseMatrix<int,ScalarType> Beta(Nc,Nc,true);
1508 
1509  computeRitzIndex();
1510 
1511  for( int i=0; i<Nc; ++i )
1512  {
1513  if( d_ritzIndex[i] == 0 )
1514  {
1515  Alpha(i,i) = d_alphar[i];
1516  Beta(i,i) = d_betar[i];
1517  }
1518  else if( d_ritzIndex[i] == 1 )
1519  {
1520  Alpha(i,i) = d_alphar[i];
1521  Alpha(i,i+1) = d_alphai[i];
1522  Beta(i,i) = d_betar[i];
1523  }
1524  else
1525  {
1526  Alpha(i,i-1) = d_alphai[i];
1527  Alpha(i,i) = d_alphar[i];
1528  Beta(i,i) = d_betar[i];
1529  }
1530  }
1531 
1532  int err;
1533 
1534  // Multiply the eigenvectors by alpha
1535  err = X_alpha.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, ST::one(), X, Alpha, ST::zero() );
1536  std::stringstream astream;
1537  astream << "GeneralizedDavidson::ScaleEigenvectors: multiply returned error code " << err;
1538  TEUCHOS_TEST_FOR_EXCEPTION( err!=0, std::runtime_error, astream.str() );
1539 
1540  // Multiply the eigenvectors by beta
1541  err = X_beta.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, ST::one(), X, Beta, ST::zero() );
1542  std::stringstream bstream;
1543  bstream << "GeneralizedDavidson::ScaleEigenvectors: multiply returned error code " << err;
1544  TEUCHOS_TEST_FOR_EXCEPTION( err!=0, std::runtime_error, bstream.str() );
1545 }
1546 
1547 //---------------------------------------------------------------------------//
1548 // Compute residual
1549 //---------------------------------------------------------------------------//
1550 template <class ScalarType, class MV, class OP>
1552 {
1553  computeRitzIndex();
1554 
1555  // Determine how many residual vectors need to be computed
1556  d_residualSize = std::max( d_blockSize, d_NEV );
1557  if( d_curDim < d_residualSize )
1558  {
1559  d_residualSize = d_curDim;
1560  }
1561  else if( d_ritzIndex[d_residualSize-1] == 1 )
1562  {
1563  d_residualSize++;
1564  }
1565 
1566  // Get indices of all valid residual vectors
1567  std::vector<int> residualIndices(d_residualSize);
1568  for( int i=0; i<d_residualSize; ++i )
1569  residualIndices[i] = i;
1570 
1571  // X will store (right) eigenvectors of projected system
1572  Teuchos::SerialDenseMatrix<int,ScalarType> X(Teuchos::Copy,*d_Z,d_curDim,d_curDim);
1573 
1574  // Get eigenvectors of projected problem -- computed from previous Schur decomposition
1575  computeProjectedEigenvectors( X );
1576 
1577  // X_alpha and X_beta will be eigenvectors right-multiplied by alpha, beta (which are quasi-diagonal portions of S,T)
1578  Teuchos::SerialDenseMatrix<int,ScalarType> X_alpha(d_curDim,d_residualSize);
1579  Teuchos::SerialDenseMatrix<int,ScalarType> X_beta(d_curDim,d_residualSize);
1580 
1581  // X_wanted is the wanted portion of X
1582  Teuchos::SerialDenseMatrix<int,ScalarType> X_wanted(Teuchos::View, X, d_curDim, d_residualSize);
1583 
1584  // Scale Eigenvectors by alpha or beta
1585  scaleEigenvectors( X_wanted, X_alpha, X_beta );
1586 
1587  // Get view of residual vector(s)
1588  RCP<MV> R_active = MVT::CloneViewNonConst( *d_R, residualIndices );
1589 
1590  // View of active portion of AV,BV
1591  std::vector<int> activeIndices(d_curDim);
1592  for( int i=0; i<d_curDim; ++i )
1593  activeIndices[i]=i;
1594 
1595  // Compute residual
1596  RCP<const MV> AV_active = MVT::CloneView( *d_AV, activeIndices );
1597  MVT::MvTimesMatAddMv(ST::one(),*AV_active, X_beta, ST::zero(),*R_active);
1598 
1599  if( d_haveB )
1600  {
1601  RCP<const MV> BV_active = MVT::CloneView( *d_BV, activeIndices );
1602  MVT::MvTimesMatAddMv(ST::one(),*BV_active, X_alpha,-ST::one(), *R_active);
1603  }
1604  else
1605  {
1606  RCP<const MV> V_active = MVT::CloneView( *d_V, activeIndices );
1607  MVT::MvTimesMatAddMv(ST::one(),*V_active, X_alpha,-ST::one(), *R_active);
1608  }
1609 
1610  /* Apply a scaling to the residual
1611  * For generalized eigenvalue problems, LAPACK scales eigenvectors
1612  * to have unit length in the infinity norm, we want them to have unit
1613  * length in the 2-norm. For conjugate pairs, the scaling is such that
1614  * |xr|^2 + |xi|^2 = 1
1615  * Additionally, the residual is currently computed as r=beta*A*x-alpha*B*x
1616  * but the "standard" residual is r=A*x-(alpha/beta)*B*x, or if we want
1617  * to scale the residual by the Ritz value then it is r=(beta/alpha)*A*x-B*x
1618  * Performing the scaling this way allows us to avoid the possibility of
1619  * diving by infinity or zero if the StatusTest were allowed to handle the
1620  * scaling.
1621  */
1622  Teuchos::LAPACK<int,ScalarType> lapack;
1623  Teuchos::BLAS<int,ScalarType> blas;
1624  std::vector<MagnitudeType> resScaling(d_residualSize);
1625  for( int icol=0; icol<d_residualSize; ++icol )
1626  {
1627  if( d_ritzIndex[icol] == 0 )
1628  {
1629  MagnitudeType Xnrm = blas.NRM2( d_curDim, X_wanted[icol], 1);
1630  MagnitudeType ABscaling = d_relativeConvergence ? d_alphar[icol] : d_betar[icol];
1631  resScaling[icol] = MT::one() / (Xnrm * ABscaling);
1632  }
1633  else if( d_ritzIndex[icol] == 1 )
1634  {
1635  MagnitudeType Xnrm1 = blas.NRM2( d_curDim, X_wanted[icol], 1 );
1636  MagnitudeType Xnrm2 = blas.NRM2( d_curDim, X_wanted[icol+1], 1 );
1637  MagnitudeType Xnrm = lapack.LAPY2(Xnrm1,Xnrm2);
1638  MagnitudeType ABscaling = d_relativeConvergence ? lapack.LAPY2(d_alphar[icol],d_alphai[icol])
1639  : d_betar[icol];
1640  resScaling[icol] = MT::one() / (Xnrm * ABscaling);
1641  resScaling[icol+1] = MT::one() / (Xnrm * ABscaling);
1642  }
1643  }
1644  MVT::MvScale( *R_active, resScaling );
1645 
1646  // Compute residual norms
1647  d_resNorms.resize(d_residualSize);
1648  MVT::MvNorm(*R_active,d_resNorms);
1649 
1650  // If Ritz value i is real, then the corresponding residual vector
1651  // is the true residual
1652  // If Ritz values i and i+1 form a conjugate pair, then the
1653  // corresponding residual vectors are the real and imaginary components
1654  // of the residual. Adjust the residual norms appropriately...
1655  for( int i=0; i<d_residualSize; ++i )
1656  {
1657  if( d_ritzIndex[i] == 1 )
1658  {
1659  MagnitudeType nrm = lapack.LAPY2(d_resNorms[i],d_resNorms[i+1]);
1660  d_resNorms[i] = nrm;
1661  d_resNorms[i+1] = nrm;
1662  }
1663  }
1664 
1665  // Evaluate with status test
1666  d_tester->checkStatus(this);
1667 
1668  // Determine which residual vectors should be used for subspace expansion
1669  if( d_useLeading || d_blockSize >= d_NEV )
1670  {
1671  d_expansionSize=d_blockSize;
1672  if( d_ritzIndex[d_blockSize-1]==1 )
1673  d_expansionSize++;
1674 
1675  d_expansionIndices.resize(d_expansionSize);
1676  for( int i=0; i<d_expansionSize; ++i )
1677  d_expansionIndices[i] = i;
1678  }
1679  else
1680  {
1681  std::vector<int> convergedVectors = d_tester->whichVecs();
1682 
1683  // Get index of first unconverged vector
1684  int startVec;
1685  for( startVec=0; startVec<d_residualSize; ++startVec )
1686  {
1687  if( std::find(convergedVectors.begin(),convergedVectors.end(),startVec)==convergedVectors.end() )
1688  break;
1689  }
1690 
1691  // Now get a contiguous block of indices starting at startVec
1692  // If this crosses the end of our residual vectors, take the final d_blockSize vectors
1693  int endVec = startVec + d_blockSize - 1;
1694  if( endVec > (d_residualSize-1) )
1695  {
1696  endVec = d_residualSize-1;
1697  startVec = d_residualSize-d_blockSize;
1698  }
1699 
1700  // Don't split conjugate pairs on either end of the range
1701  if( d_ritzIndex[startVec]==-1 )
1702  {
1703  startVec--;
1704  endVec--;
1705  }
1706 
1707  if( d_ritzIndex[endVec]==1 )
1708  endVec++;
1709 
1710  d_expansionSize = 1+endVec-startVec;
1711  d_expansionIndices.resize(d_expansionSize);
1712  for( int i=0; i<d_expansionSize; ++i )
1713  d_expansionIndices[i] = startVec+i;
1714  }
1715 }
1716 
1717 //---------------------------------------------------------------------------//
1718 // Print current status.
1719 //---------------------------------------------------------------------------//
1720 template <class ScalarType, class MV, class OP>
1722 {
1723  using std::endl;
1724 
1725  myout.setf(std::ios::scientific, std::ios::floatfield);
1726  myout.precision(6);
1727  myout <<endl;
1728  myout <<"================================================================================" << endl;
1729  myout << endl;
1730  myout <<" GeneralizedDavidson Solver Solver Status" << endl;
1731  myout << endl;
1732  myout <<"The solver is "<<(d_initialized ? "initialized." : "not initialized.") << endl;
1733  myout <<"The number of iterations performed is " << d_iteration << endl;
1734  myout <<"The number of operator applies performed is " << d_opApplies << endl;
1735  myout <<"The block size is " << d_expansionSize << endl;
1736  myout <<"The current basis size is " << d_curDim << endl;
1737  myout <<"The number of requested eigenvalues is " << d_NEV << endl;
1738  myout <<"The number of converged values is " << d_tester->howMany() << endl;
1739  myout << endl;
1740 
1741  myout.setf(std::ios_base::right, std::ios_base::adjustfield);
1742 
1743  if( d_initialized )
1744  {
1745  myout << "CURRENT RITZ VALUES" << endl;
1746 
1747  myout << std::setw(24) << "Ritz Value"
1748  << std::setw(30) << "Residual Norm" << endl;
1749  myout << "--------------------------------------------------------------------------------" << endl;
1750  if( d_residualSize > 0 )
1751  {
1752  std::vector<MagnitudeType> realRitz(d_curDim), imagRitz(d_curDim);
1753  std::vector< Value<ScalarType> > ritzVals = getRitzValues();
1754  for( int i=0; i<d_curDim; ++i )
1755  {
1756  realRitz[i] = ritzVals[i].realpart;
1757  imagRitz[i] = ritzVals[i].imagpart;
1758  }
1759  std::vector<int> permvec;
1760  sortValues( realRitz, imagRitz, permvec, d_curDim );
1761 
1762  int numToPrint = std::max( d_numToPrint, d_NEV );
1763  numToPrint = std::min( d_curDim, numToPrint );
1764 
1765  // Because the sort manager does not use a stable sort, occasionally
1766  // the portion of a conjugate pair with negative imaginary part will be placed
1767  // first...in that case the following will not give the usual expected behavior
1768  // and an extra value will be printed. This is only an issue with the output
1769  // format because the actually sorting of Schur forms is guaranteed to be stable.
1770  if( d_ritzIndex[permvec[numToPrint-1]] != 0 )
1771  numToPrint++;
1772 
1773  int i=0;
1774  while( i<numToPrint )
1775  {
1776  if( imagRitz[i] == ST::zero() )
1777  {
1778  myout << std::setw(15) << realRitz[i];
1779  myout << " + i" << std::setw(15) << ST::magnitude( imagRitz[i] );
1780  if( i < d_residualSize )
1781  myout << std::setw(20) << d_resNorms[permvec[i]] << endl;
1782  else
1783  myout << " Not Computed" << endl;
1784 
1785  i++;
1786  }
1787  else
1788  {
1789  // Positive imaginary part
1790  myout << std::setw(15) << realRitz[i];
1791  myout << " + i" << std::setw(15) << ST::magnitude( imagRitz[i] );
1792  if( i < d_residualSize )
1793  myout << std::setw(20) << d_resNorms[permvec[i]] << endl;
1794  else
1795  myout << " Not Computed" << endl;
1796 
1797  // Negative imaginary part
1798  myout << std::setw(15) << realRitz[i];
1799  myout << " - i" << std::setw(15) << ST::magnitude( imagRitz[i] );
1800  if( i < d_residualSize )
1801  myout << std::setw(20) << d_resNorms[permvec[i]] << endl;
1802  else
1803  myout << " Not Computed" << endl;
1804 
1805  i+=2;
1806  }
1807  }
1808  }
1809  else
1810  {
1811  myout << std::setw(20) << "[ NONE COMPUTED ]" << endl;
1812  }
1813  }
1814  myout << endl;
1815  myout << "================================================================================" << endl;
1816  myout << endl;
1817 }
1818 
1819 } // namespace Anasazi
1820 
1821 #endif // ANASAZI_GENERALIZED_DAVIDSON_HPP
1822 
void setStatusTest(RCP< StatusTest< ScalarType, MV, OP > > tester)
Set status test.
RCP< MV > V
Orthonormal basis for search subspace.
std::vector< MagnitudeType > getRitzRes2Norms()
Get the current Ritz residual norms (2-norm)
RCP< StatusTest< ScalarType, MV, OP > > getStatusTest() const
Get status test.
Teuchos::ScalarTraits< ScalarType >::magnitudeType imagpart
The imaginary component of the eigenvalue.
This class defines the interface required by an eigensolver and status test class to compute solution...
int getBlockSize() const
Get block size.
Solves eigenvalue problem using generalized Davidson method.
std::vector< Value< ScalarType > > getRitzValues()
Get the current Ritz values.
RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > T
Right quasi upper triangular matrix from QZ decomposition of (VAV,VBV)
Virtual base class which defines basic traits for the operator type.
GeneralizedDavidsonState< ScalarType, MV > getState()
Get the current state of the eigensolver.
RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > VAV
Projection of A onto V.
void setBlockSize(int blockSize)
Set block size.
int curDim
The current subspace dimension.
std::vector< int > getRitzIndex()
Get the current Ritz index vector.
const Eigenproblem< ScalarType, MV, OP > & getProblem() const
Get eigenproblem.
Pure virtual base class which describes the basic interface to the iterative eigensolver.
Namespace Anasazi contains the classes, structs, enums and utilities used by the Anasazi package...
void iterate()
Solves the eigenvalue problem.
std::vector< MagnitudeType > getRes2Norms()
Get the current residual norms (2-norm)
Structure to contain pointers to GeneralizedDavidson state variables.
Abstract class definition for Anasazi Output Managers.
std::vector< Value< ScalarType > > eVals
Vector of generalized eigenvalues.
void setSize(int blockSize, int maxSubDim)
Set problem size.
Abstract base class which defines the interface required by an eigensolver and status test class to c...
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 auxilliary vectors.
std::vector< MagnitudeType > getResNorms()
Get the current residual norms (w.r.t. norm defined by OrthoManager)
Templated virtual class for providing orthogonalization/orthonormalization methods.
RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > Q
Left generalized Schur vectors from QZ decomposition of (VAV,VBV)
Traits class which defines basic operations on multivectors.
This struct is used for storing eigenvalues and Ritz values, as a pair of real values.
Anasazi header file which uses auto-configuration information to include necessary C++ headers...
Virtual base class which defines the interface between an eigensolver and a class whose job is the so...
bool isInitialized() const
Query if the solver is in an initialized state.
GeneralizedDavidson(const RCP< Eigenproblem< ScalarType, MV, OP > > &problem, const RCP< SortManager< MagnitudeType > > &sortman, const RCP< OutputManager< ScalarType > > &outputman, const RCP< StatusTest< ScalarType, MV, OP > > &tester, const RCP< OrthoManager< ScalarType, MV > > &orthoman, Teuchos::ParameterList &pl)
Constructor.
int getCurSubspaceDim() const
Get current subspace dimension.
RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > S
Left quasi upper triangular matrix from QZ decomposition of (VAV,VBV)
RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > VBV
Projection of B onto V.
void initialize()
Initialize the eigenvalue problem.
int getNumIters() const
Get number of iterations.
RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > Z
Right generalized Schur vectors from QZ decomposition of (VAV,VBV)
void setAuxVecs(const Teuchos::Array< RCP< const MV > > &auxVecs)
Set auxilliary vectors.
Types and exceptions used within Anasazi solvers and interfaces.
Anasazi&#39;s templated virtual class for providing routines for orthogonalization and orthonormalization...
Anasazi&#39;s templated pure virtual class for managing the sorting of approximate eigenvalues computed b...
RCP< const MV > getRitzVectors()
Get the current Ritz vectors.
void currentStatus(std::ostream &myout)
Print current status of solver.
Teuchos::ScalarTraits< ScalarType >::magnitudeType realpart
The real component of the eigenvalue.
Common interface of stopping criteria for Anasazi&#39;s solvers.
std::vector< int > getBlockIndex() const
Get indices of current block.
The Eigensolver is a templated virtual base class that defines the basic interface that any eigensolv...
int getMaxSubspaceDim() const
Get maximum subspace dimension.
Declaration and definition of Anasazi::StatusTest.
void resetNumIters()
Reset the number of iterations.