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

The Eigensolver is a templated virtual base class that defines the basic interface that any eigensolver will support. More...

#include <AnasaziEigensolverDecl.hpp>

Inheritance diagram for Anasazi::Eigensolver< ScalarType, MV, OP >:
Anasazi::BlockDavidson< ScalarType, MV, OP > Anasazi::BlockKrylovSchur< ScalarType, MV, OP > Anasazi::Experimental::TraceMinBase< ScalarType, MV, OP > Anasazi::GeneralizedDavidson< ScalarType, MV, OP > Anasazi::LOBPCG< ScalarType, MV, OP > Anasazi::RTRBase< ScalarType, MV, OP > Anasazi::Experimental::TraceMin< ScalarType, MV, OP > Anasazi::Experimental::TraceMinDavidson< ScalarType, MV, OP > Anasazi::IRTR< ScalarType, MV, OP > Anasazi::SIRTR< ScalarType, MV, OP >

Public Member Functions

Constructors/Destructor
 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...
 
Solver methods
virtual void iterate ()=0
 This method performs eigensolvers iterations until the status test indicates the need to stop or an error occurs (in which case, an exception is thrown). More...
 
virtual void initialize ()=0
 Initialize the solver with the initial vectors from the eigenproblem or random data. More...
 
Status methods
virtual int getNumIters () const =0
 Get the current iteration count. More...
 
virtual void resetNumIters ()=0
 Reset the iteration count. More...
 
virtual Teuchos::RCP< const MV > getRitzVectors ()=0
 Get the Ritz vectors from the previous iteration. These are indexed using getRitzIndex(). More...
 
virtual std::vector< Value
< ScalarType > > 
getRitzValues ()=0
 Get the Ritz values from the previous iteration. More...
 
virtual std::vector< int > getRitzIndex ()=0
 Get the index used for indexing the compressed storage used for Ritz vectors for real, non-Hermitian problems. More...
 
virtual std::vector< typename
Teuchos::ScalarTraits
< ScalarType >::magnitudeType > 
getResNorms ()=0
 Get the current residual norms. More...
 
virtual std::vector< typename
Teuchos::ScalarTraits
< ScalarType >::magnitudeType > 
getRes2Norms ()=0
 
virtual std::vector< typename
Teuchos::ScalarTraits
< ScalarType >::magnitudeType > 
getRitzRes2Norms ()=0
 
virtual int getCurSubspaceDim () const =0
 Get the dimension of the search subspace used to generate the current eigenvectors and eigenvalues. More...
 
virtual int getMaxSubspaceDim () const =0
 Get the maximum dimension allocated for the search subspace. More...
 
Accessor methods
virtual void setStatusTest (Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test)=0
 Set a new StatusTest for the solver. More...
 
virtual Teuchos::RCP
< StatusTest< ScalarType, MV,
OP > > 
getStatusTest () const =0
 Get the current StatusTest used by the solver. More...
 
virtual const Eigenproblem
< ScalarType, MV, OP > & 
getProblem () const =0
 Get a constant reference to the eigenvalue problem. More...
 
virtual int getBlockSize () const =0
 Get the blocksize to be used by the iterative solver in solving this eigenproblem. More...
 
virtual void setBlockSize (int blockSize)=0
 Set the blocksize to be used by the iterative solver in solving this eigenproblem. More...
 
virtual void setAuxVecs (const Teuchos::Array< Teuchos::RCP< const MV > > &auxvecs)=0
 Set the auxiliary vectors for the solver. More...
 
virtual Teuchos::Array
< Teuchos::RCP< const MV > > 
getAuxVecs () const =0
 Get the auxiliary vectors for the solver. More...
 
virtual bool isInitialized () const =0
 States whether the solver has been initialized or not. More...
 
Output methods
virtual void currentStatus (std::ostream &os)=0
 This method requests that the solver print out its current status to screen. More...
 

Detailed Description

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

The Eigensolver is a templated virtual base class that defines the basic interface that any eigensolver will support.

This interface is mainly concerned with providing a set of eigensolver status method that can be requested from any eigensolver by an StatusTest object.

Examples:
LOBPCGCustomStatusTest.cpp.

Definition at line 67 of file AnasaziEigensolver.hpp.

Constructor & Destructor Documentation

template<class ScalarType, class MV, class OP>
Anasazi::Eigensolver< ScalarType, MV, OP >::Eigensolver ( )
inline

Default Constructor.

Definition at line 75 of file AnasaziEigensolver.hpp.

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

This constructor, implemented by all Anasazi eigensolvers, takes an Anasazi::Eigenproblem, Anasazi::SortManager, Anasazi::OutputManager, and Teuchos::ParameterList as input. These four arguments are sufficient enough for constructing any Anasazi::Eigensolver object.

template<class ScalarType, class MV, class OP>
virtual Anasazi::Eigensolver< ScalarType, MV, OP >::~Eigensolver ( )
inlinevirtual

Destructor.

Definition at line 90 of file AnasaziEigensolver.hpp.

Member Function Documentation

template<class ScalarType, class MV, class OP>
virtual void Anasazi::Eigensolver< ScalarType, MV, OP >::iterate ( )
pure virtual
template<class ScalarType, class MV, class OP>
virtual void Anasazi::Eigensolver< ScalarType, MV, OP >::initialize ( )
pure virtual
template<class ScalarType, class MV, class OP>
virtual int Anasazi::Eigensolver< ScalarType, MV, OP >::getNumIters ( ) const
pure virtual
template<class ScalarType, class MV, class OP>
virtual void Anasazi::Eigensolver< ScalarType, MV, OP >::resetNumIters ( )
pure virtual
template<class ScalarType, class MV, class OP>
virtual Teuchos::RCP<const MV> Anasazi::Eigensolver< ScalarType, MV, OP >::getRitzVectors ( )
pure virtual
template<class ScalarType, class MV, class OP>
virtual std::vector<Value<ScalarType> > Anasazi::Eigensolver< ScalarType, MV, OP >::getRitzValues ( )
pure virtual
template<class ScalarType, class MV, class OP>
virtual std::vector<int> Anasazi::Eigensolver< ScalarType, MV, OP >::getRitzIndex ( )
pure virtual

Get the index used for indexing the compressed storage used for Ritz vectors for real, non-Hermitian problems.

index has length numVecs, where each entry is 0, +1, or -1. These have the following interpretation:

  • index[i] == 0: signifies that the corresponding eigenvector is stored as the i column of Evecs. This will usually be the case when ScalarType is complex, an eigenproblem is Hermitian, or a real, non-Hermitian eigenproblem has a real eigenvector.
  • index[i] == +1: signifies that the corresponding eigenvector is stored in two vectors: the real part in the i column of Evecs and the positive imaginary part in the i+1 column of Evecs.
  • index[i] == -1: signifies that the corresponding eigenvector is stored in two vectors: the real part in the i-1 column of Evecs and the negative imaginary part in the i column of Evecs

Implemented in Anasazi::LOBPCG< ScalarType, MV, OP >, Anasazi::Experimental::TraceMinBase< ScalarType, MV, OP >, Anasazi::RTRBase< ScalarType, MV, OP >, Anasazi::BlockDavidson< ScalarType, MV, OP >, Anasazi::BlockKrylovSchur< ScalarType, MV, OP >, and Anasazi::GeneralizedDavidson< ScalarType, MV, OP >.

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

Get the current residual norms.

Returns
A vector of length blockSize containing the norms of the residuals, according to the orthogonalization manager norm() method.

Implemented in Anasazi::LOBPCG< ScalarType, MV, OP >, Anasazi::Experimental::TraceMinBase< ScalarType, MV, OP >, Anasazi::RTRBase< ScalarType, MV, OP >, Anasazi::BlockDavidson< ScalarType, MV, OP >, Anasazi::BlockKrylovSchur< ScalarType, MV, OP >, and Anasazi::GeneralizedDavidson< ScalarType, MV, OP >.

template<class ScalarType, class MV, class OP>
virtual std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> Anasazi::Eigensolver< ScalarType, MV, OP >::getRes2Norms ( )
pure virtual
template<class ScalarType, class MV, class OP>
virtual std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> Anasazi::Eigensolver< ScalarType, MV, OP >::getRitzRes2Norms ( )
pure virtual
template<class ScalarType, class MV, class OP>
virtual int Anasazi::Eigensolver< ScalarType, MV, OP >::getCurSubspaceDim ( ) const
pure virtual
template<class ScalarType, class MV, class OP>
virtual int Anasazi::Eigensolver< ScalarType, MV, OP >::getMaxSubspaceDim ( ) const
pure virtual
template<class ScalarType, class MV, class OP>
virtual void Anasazi::Eigensolver< ScalarType, MV, OP >::setStatusTest ( Teuchos::RCP< StatusTest< ScalarType, MV, OP > >  test)
pure virtual
template<class ScalarType, class MV, class OP>
virtual Teuchos::RCP<StatusTest<ScalarType,MV,OP> > Anasazi::Eigensolver< ScalarType, MV, OP >::getStatusTest ( ) const
pure virtual
template<class ScalarType, class MV, class OP>
virtual const Eigenproblem<ScalarType,MV,OP>& Anasazi::Eigensolver< ScalarType, MV, OP >::getProblem ( ) const
pure virtual
template<class ScalarType, class MV, class OP>
virtual int Anasazi::Eigensolver< ScalarType, MV, OP >::getBlockSize ( ) const
pure virtual
template<class ScalarType, class MV, class OP>
virtual void Anasazi::Eigensolver< ScalarType, MV, OP >::setBlockSize ( int  blockSize)
pure virtual
template<class ScalarType, class MV, class OP>
virtual void Anasazi::Eigensolver< ScalarType, MV, OP >::setAuxVecs ( const Teuchos::Array< Teuchos::RCP< const MV > > &  auxvecs)
pure virtual
template<class ScalarType, class MV, class OP>
virtual Teuchos::Array<Teuchos::RCP<const MV> > Anasazi::Eigensolver< ScalarType, MV, OP >::getAuxVecs ( ) const
pure virtual
template<class ScalarType, class MV, class OP>
virtual bool Anasazi::Eigensolver< ScalarType, MV, OP >::isInitialized ( ) const
pure virtual
template<class ScalarType, class MV, class OP>
virtual void Anasazi::Eigensolver< ScalarType, MV, OP >::currentStatus ( std::ostream &  os)
pure virtual

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