Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | List of all members
Anasazi::GeneralizedDavidson< ScalarType, MV, OP > Class Template Reference

Solves eigenvalue problem using generalized Davidson method. More...

#include <AnasaziGeneralizedDavidson.hpp>

Inheritance diagram for Anasazi::GeneralizedDavidson< ScalarType, MV, OP >:
Anasazi::Eigensolver< ScalarType, MV, OP >

Public Member Functions

 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. More...
 
void iterate ()
 Solves the eigenvalue problem. More...
 
void initialize ()
 Initialize the eigenvalue problem. More...
 
void initialize (GeneralizedDavidsonState< ScalarType, MV > &state)
 Initialize solver from state. More...
 
int getNumIters () const
 Get number of iterations. More...
 
void resetNumIters ()
 Reset the number of iterations. More...
 
RCP< const MV > getRitzVectors ()
 Get the current Ritz vectors. More...
 
std::vector< Value< ScalarType > > getRitzValues ()
 Get the current Ritz values. More...
 
std::vector< int > getRitzIndex ()
 Get the current Ritz index vector. More...
 
std::vector< int > getBlockIndex () const
 Get indices of current block. More...
 
std::vector< MagnitudeType > getResNorms ()
 Get the current residual norms (w.r.t. norm defined by OrthoManager) More...
 
std::vector< MagnitudeType > getResNorms (int numWanted)
 Get the current residual norms (w.r.t. norm defined by OrthoManager) More...
 
std::vector< MagnitudeType > getRes2Norms ()
 Get the current residual norms (2-norm) More...
 
std::vector< MagnitudeType > getRitzRes2Norms ()
 Get the current Ritz residual norms (2-norm) More...
 
int getCurSubspaceDim () const
 Get current subspace dimension. More...
 
int getMaxSubspaceDim () const
 Get maximum subspace dimension. More...
 
void setStatusTest (RCP< StatusTest< ScalarType, MV, OP > > tester)
 Set status test. More...
 
RCP< StatusTest< ScalarType,
MV, OP > > 
getStatusTest () const
 Get status test. More...
 
const Eigenproblem< ScalarType,
MV, OP > & 
getProblem () const
 Get eigenproblem. More...
 
int getBlockSize () const
 Get block size. More...
 
void setBlockSize (int blockSize)
 Set block size. More...
 
void setSize (int blockSize, int maxSubDim)
 Set problem size. More...
 
Teuchos::Array< RCP< const MV > > getAuxVecs () const
 Get the auxilliary vectors. More...
 
void setAuxVecs (const Teuchos::Array< RCP< const MV > > &auxVecs)
 Set auxilliary vectors. More...
 
bool isInitialized () const
 Query if the solver is in an initialized state. More...
 
void currentStatus (std::ostream &myout)
 Print current status of solver. More...
 
GeneralizedDavidsonState
< ScalarType, MV > 
getState ()
 Get the current state of the eigensolver. More...
 
void sortProblem (int numWanted)
 
- Public Member Functions inherited from Anasazi::Eigensolver< ScalarType, MV, OP >
 Eigensolver ()
 Default Constructor. More...
 
 Eigensolver (const Teuchos::RCP< Eigenproblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< SortManager< ScalarType > > &sorter, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, const Teuchos::RCP< OrthoManager< ScalarType, MV > > &ortho, Teuchos::ParameterList &params)
 Basic Constructor. More...
 
virtual ~Eigensolver ()
 Destructor. More...
 

Detailed Description

template<class ScalarType, class MV, class OP>
class Anasazi::GeneralizedDavidson< ScalarType, MV, OP >

Solves eigenvalue problem using generalized Davidson method.

This class searches for a few eigenvalues and corresponding eigenvectors for either a standard eigenvalue problem $Ax=\lambda x$ or a generalized eigenvalue problem $Ax=\lambda B x$ Note that unlike some other solvers, the generalized Davidson method places no restrictions on either matrix in a generalized eigenvalue problem.

Tips for preconditioning: A good preconditioner usually approximates $(A-\sigma I)^{-1}$ or $(A-\sigma B)^{-1}$, where $\sigma$ is close to the target eigenvalue. When searching for largest magnitude eigenvalues, selecting a preconditioner $P^{-1} \approx B^{-1}$ usually works well and when searching for smallest magnitude eigenvalues selecting $P^{-1} \approx A^{-1}$ is usually appropriate.

This class is currently only implemented for real scalar types (i.e. float, double).

Definition at line 141 of file AnasaziGeneralizedDavidson.hpp.

Constructor & Destructor Documentation

template<class ScalarType , class MV , class OP >
Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::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.

GeneralizedDavidson constructor with eigenproblem, parameters, and solver utilities.

Behavior of the solver is controlled by the following ParameterList entries:

  • "Block Size" – block size used by algorithm. Default: 1.
  • "Maximum Subspace Dimension" – maximum number of basis vectors for subspace. Two for standard eigenvalue problem) or three (for generalized eigenvalue problem) sets of basis vectors of this size will be required. Default: 3*problem->getNEV()*"Block Size"
  • "Initial Guess" – how should initial vector be selected: "Random" or "User". If "User," the value in problem->getInitVec() will be used. Default: "Random".
  • "Print Number of Ritz Values" – an int specifying how many Ritz values should be printed at each iteration. Default: "NEV".
  • "Relative Convergence Tolerance" – should residual be scaled by corresponding Ritz value to measure convergence. Default: "false"

Definition at line 489 of file AnasaziGeneralizedDavidson.hpp.

Member Function Documentation

template<class ScalarType , class MV , class OP >
void Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::iterate ( )
virtual

Solves the eigenvalue problem.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 567 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::initialize ( )
virtual

Initialize the eigenvalue problem.

Anything on the state that is not null is assumed to be valid. Anything not present on the state will be generated. Very limited error checking can be performed to ensure the validity of state components (e.g. we cannot verify that state.AV actually corresponds to A*state.V), so this function should be used carefully.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 887 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::initialize ( GeneralizedDavidsonState< ScalarType, MV > &  state)

Initialize solver from state.

Definition at line 720 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
int Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getNumIters ( ) const
inlinevirtual

Get number of iterations.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 204 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::resetNumIters ( )
inlinevirtual

Reset the number of iterations.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 209 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
RCP<const MV> Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getRitzVectors ( )
inlinevirtual

Get the current Ritz vectors.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 214 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
std::vector< Value< ScalarType > > Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getRitzValues ( )
virtual

Get the current Ritz values.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 926 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
std::vector<int> Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getRitzIndex ( )
inlinevirtual

Get the current Ritz index vector.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 229 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
std::vector<int> Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getBlockIndex ( ) const
inline

Get indices of current block.

Number of entries is equal to getBlockSize()

Definition at line 241 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getResNorms ( )
virtual

Get the current residual norms (w.r.t. norm defined by OrthoManager)

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 898 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getResNorms ( int  numWanted)

Get the current residual norms (w.r.t. norm defined by OrthoManager)

Definition at line 908 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
std::vector<MagnitudeType> Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getRes2Norms ( )
inlinevirtual

Get the current residual norms (2-norm)

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 259 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
std::vector<MagnitudeType> Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getRitzRes2Norms ( )
inlinevirtual

Get the current Ritz residual norms (2-norm)

GeneralizedDavidson doesn't compute Ritz residual norms so this is equivalent to calling getRes2Norms()

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 267 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
int Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getCurSubspaceDim ( ) const
inlinevirtual

Get current subspace dimension.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 272 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
int Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getMaxSubspaceDim ( ) const
inlinevirtual

Get maximum subspace dimension.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 277 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::setStatusTest ( RCP< StatusTest< ScalarType, MV, OP > >  tester)
inlinevirtual

Set status test.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 282 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
RCP<StatusTest<ScalarType,MV,OP> > Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getStatusTest ( ) const
inlinevirtual

Get status test.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 287 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
const Eigenproblem<ScalarType,MV,OP>& Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getProblem ( ) const
inlinevirtual

Get eigenproblem.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 292 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
int Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getBlockSize ( ) const
inlinevirtual

Get block size.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 297 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::setBlockSize ( int  blockSize)
virtual

Set block size.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 645 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::setSize ( int  blockSize,
int  maxSubDim 
)

Set problem size.

Definition at line 654 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::Array< RCP<const MV> > Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getAuxVecs ( ) const
inlinevirtual

Get the auxilliary vectors.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 312 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::setAuxVecs ( const Teuchos::Array< RCP< const MV > > &  auxVecs)
virtual

Set auxilliary vectors.

Manually setting the auxilliary vectors invalidates the current state of the solver. Reuse of any components of the solver requires extracting the state, orthogonalizing V against the aux vecs and reinitializing.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 954 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
bool Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::isInitialized ( ) const
inlinevirtual

Query if the solver is in an initialized state.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 326 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::currentStatus ( std::ostream &  myout)
virtual

Print current status of solver.

Implements Anasazi::Eigensolver< ScalarType, MV, OP >.

Definition at line 1722 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
GeneralizedDavidsonState< ScalarType, MV > Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::getState ( )

Get the current state of the eigensolver.

Definition at line 624 of file AnasaziGeneralizedDavidson.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::GeneralizedDavidson< ScalarType, MV, OP >::sortProblem ( int  numWanted)

Reorder Schur form, bringing wanted values to front

Definition at line 974 of file AnasaziGeneralizedDavidson.hpp.


The documentation for this class was generated from the following file: