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

This class is an abstract base class for Implicit Riemannian Trust-Region based eigensolvers. The class provides the interfaces shared by the IRTR solvers (e.g., getState() and initialize()) as well as the shared implementations (e.g., inner products). More...

#include <AnasaziRTRBase.hpp>

Inheritance diagram for Anasazi::RTRBase< ScalarType, MV, OP >:
Anasazi::Eigensolver< ScalarType, MV, OP > Anasazi::IRTR< ScalarType, MV, OP > Anasazi::SIRTR< ScalarType, MV, OP >

Public Member Functions

Constructor/Destructor
 RTRBase (const Teuchos::RCP< Eigenproblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< SortManager< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > > &sorter, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, const Teuchos::RCP< GenOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList &params, const std::string &solverLabel, bool skinnySolver)
 RTRBase constructor with eigenproblem, solver utilities, and parameter list of solver options. More...
 
virtual ~RTRBase ()
 RTRBase destructor More...
 
Solver methods
virtual void iterate ()=0
 This method performs RTR iterations until the status test indicates the need to stop or an error occurs (in which case, an exception is thrown). More...
 
void initialize (RTRState< ScalarType, MV > &newstate)
 Initialize the solver to an iterate, optionally providing the Ritz values and residual. More...
 
void initialize ()
 Initialize the solver with the initial vectors from the eigenproblem or random data. More...
 
bool isInitialized () const
 Indicates whether the solver has been initialized or not. More...
 
RTRState< ScalarType, MV > getState () const
 Get the current state of the eigensolver. More...
 
Status methods
int getNumIters () const
 Get the current iteration count. More...
 
void resetNumIters ()
 Reset the iteration count. More...
 
Teuchos::RCP< const MV > getRitzVectors ()
 Get the Ritz vectors from the previous iteration. More...
 
std::vector< Value< ScalarType > > getRitzValues ()
 Get the Ritz values from the previous iteration. More...
 
std::vector< int > getRitzIndex ()
 Get the index used for extracting Ritz vectors from getRitzVectors(). More...
 
std::vector< typename
Teuchos::ScalarTraits
< ScalarType >::magnitudeType > 
getResNorms ()
 Get the current residual norms. More...
 
std::vector< typename
Teuchos::ScalarTraits
< ScalarType >::magnitudeType > 
getRes2Norms ()
 Get the current residual 2-norms. More...
 
std::vector< typename
Teuchos::ScalarTraits
< ScalarType >::magnitudeType > 
getRitzRes2Norms ()
 Get the 2-norms of the Ritz residuals. More...
 
int getCurSubspaceDim () const
 Get the dimension of the search subspace used to generate the current eigenvectors and eigenvalues. More...
 
int getMaxSubspaceDim () const
 Get the maximum dimension allocated for the search subspace. For RTR, this always returns getBlockSize(). More...
 
Accessor routines from Eigensolver
void setStatusTest (Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test)
 Set a new StatusTest for the solver. More...
 
Teuchos::RCP< StatusTest
< ScalarType, MV, OP > > 
getStatusTest () const
 Get the current StatusTest used by the solver. More...
 
const Eigenproblem< ScalarType,
MV, OP > & 
getProblem () const
 Get a constant reference to the eigenvalue problem. More...
 
void setBlockSize (int blockSize)
 Set the blocksize to be used by the iterative solver in solving this eigenproblem. More...
 
int getBlockSize () const
 Get the blocksize to be used by the iterative solver in solving this eigenproblem. More...
 
void setAuxVecs (const Teuchos::Array< Teuchos::RCP< const MV > > &auxvecs)
 Set the auxiliary vectors for the solver. More...
 
Teuchos::Array< Teuchos::RCP
< const MV > > 
getAuxVecs () const
 Get the current auxiliary vectors. More...
 
Output methods
virtual void currentStatus (std::ostream &os)
 This method requests that the solver print out its current status to screen. More...
 
- 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::RTRBase< ScalarType, MV, OP >

This class is an abstract base class for Implicit Riemannian Trust-Region based eigensolvers. The class provides the interfaces shared by the IRTR solvers (e.g., getState() and initialize()) as well as the shared implementations (e.g., inner products).

IRTR eigensolvers are capable of solving symmetric/Hermitian eigenvalue problems. These solvers may be used to compute either the leftmost (smallest real, "SR") or rightmost (largest real, "LR") eigenvalues. For more information, see the publications at the RTR eigensolvers page.

This class is abstract and objects cannot be instantiated. Instead, instantiate one of the concrete derived classes: IRTR and SIRTR, the caching and non-caching implementations of this solver. The main difference between these solver is the memory allocated by the solvers in support of the IRTR iteration.

The reduction in memory usage is effected by eliminating the caching of operator applications. This also results in a reduction in vector arithmetic required to maintain these caches. The cost is an increase in the number of operator applications. For inexpensive operator applications, SIRTR should provide better performance over IRTR. As the operator applications becomes more expensive, the performance scale tips towards the IRTR solver. Note, the trajectory of both solvers is identical in exact arithmetic. However, the effects of round-off error in the cached results mean that some difference between the solvers may exist. This effect is seen when a large number of iterations are required to solve the trust-region subproblem in solveTRSubproblem(). Also note, the inclusion of auxiliary vectors increases the memory requirements of these solvers linearly with the number of auxiliary vectors. The required storage is listed in the following table:

Number of vectors (bS == blockSize())
SolverBase requirementGeneralized/B != nullPreconditionedGeneralized and Preconditioned
IRTR10*bS11*bS12*bS13*bS
SIRTR6*bS7*bS7*bS8*bS
Author
Chris Baker

Definition at line 202 of file AnasaziRTRBase.hpp.

Constructor & Destructor Documentation

template<class ScalarType , class MV , class OP >
Anasazi::RTRBase< ScalarType, MV, OP >::RTRBase ( const Teuchos::RCP< Eigenproblem< ScalarType, MV, OP > > &  problem,
const Teuchos::RCP< SortManager< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > > &  sorter,
const Teuchos::RCP< OutputManager< ScalarType > > &  printer,
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &  tester,
const Teuchos::RCP< GenOrthoManager< ScalarType, MV, OP > > &  ortho,
Teuchos::ParameterList params,
const std::string &  solverLabel,
bool  skinnySolver 
)

RTRBase constructor with eigenproblem, solver utilities, and parameter list of solver options.

The RTRBase class is abstract and cannot be instantiated; this constructor is called by derived classes IRTR and RTR.

Definition at line 608 of file AnasaziRTRBase.hpp.

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

RTRBase destructor

Definition at line 223 of file AnasaziRTRBase.hpp.

Member Function Documentation

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

This method performs RTR iterations until the status test indicates the need to stop or an error occurs (in which case, an exception is thrown).

iterate() will first determine whether the solver is initialized; if not, it will call initialize() using default arguments. After initialization, the solver performs RTR iterations until the status test evaluates as Passed, at which point the method returns to the caller.

The RTR iteration proceeds as follows:

  1. the trust-region subproblem at X is solved for update Eta via a call to solveTRSubproblem()
  2. the new iterate is the Ritz vectors with respect to X+Eta
  3. the eigenproblem residuals are formed with respect to the new iterate

The status test is queried at the beginning of the iteration.

Possible exceptions thrown include std::logic_error, std::invalid_argument or one of the RTR-specific exceptions.

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

Implemented in Anasazi::IRTR< ScalarType, MV, OP >, and Anasazi::SIRTR< ScalarType, MV, OP >.

template<class ScalarType , class MV , class OP >
void Anasazi::RTRBase< ScalarType, MV, OP >::initialize ( RTRState< ScalarType, MV > &  newstate)

Initialize the solver to an iterate, optionally providing the Ritz values and residual.

The RTR eigensolver contains a certain amount of state relating to the current iterate.

initialize() gives the user the opportunity to manually set these, although this must be done with caution, abiding by the rules given below. All notions of orthogonality and orthonormality are derived from the inner product specified by the orthogonalization manager.

Postcondition

The user has the option of specifying any component of the state using initialize(). However, these arguments are assumed to match the post-conditions specified under isInitialized(). Any component of the state (i.e., AX) not given to initialize() will be generated.

If the Ritz values relative to newstate.X are passed in newstate.T, then newstate.X is assume to contain Ritz vectors, i.e., newstate.T must be B-orthonormal and it must partially diagonalize A.

Definition at line 1071 of file AnasaziRTRBase.hpp.

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

Initialize the solver with the initial vectors from the eigenproblem or random data.

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

Definition at line 1435 of file AnasaziRTRBase.hpp.

template<class ScalarType , class MV , class OP >
bool Anasazi::RTRBase< ScalarType, MV, OP >::isInitialized ( ) const
virtual

Indicates whether the solver has been initialized or not.

Returns
bool indicating the state of the solver.
Postcondition
If isInitialized() == true:
  • X is orthogonal to auxiliary vectors and has orthonormal columns
  • AX == A*X
  • BX == B*X if B != Teuchos::null
    Otherwise, BX == Teuchos::null
  • getRitzValues() returns the sorted Ritz values with respect to X
  • getResidualVecs() returns the residual vectors with respect to X

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

Definition at line 1781 of file AnasaziRTRBase.hpp.

template<class ScalarType , class MV , class OP >
RTRState< ScalarType, MV > Anasazi::RTRBase< ScalarType, MV, OP >::getState ( ) const

Get the current state of the eigensolver.

The data is only valid if isInitialized() == true.

Returns
An RTRState object containing const pointers to the current solver state.

Definition at line 1763 of file AnasaziRTRBase.hpp.

template<class ScalarType , class MV , class OP >
int Anasazi::RTRBase< ScalarType, MV, OP >::getNumIters ( ) const
virtual

Get the current iteration count.

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

Definition at line 1757 of file AnasaziRTRBase.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::RTRBase< ScalarType, MV, OP >::resetNumIters ( )
virtual

Reset the iteration count.

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

Definition at line 1751 of file AnasaziRTRBase.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP< const MV > Anasazi::RTRBase< ScalarType, MV, OP >::getRitzVectors ( )
virtual

Get the Ritz vectors from the previous iteration.

Returns
A multivector with getBlockSize() vectors containing the sorted Ritz vectors corresponding to the most significant Ritz values. The i-th vector of the return corresponds to the i-th Ritz vector; there is no need to use getRitzIndex().

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

Definition at line 1745 of file AnasaziRTRBase.hpp.

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

Get the Ritz values from the previous iteration.

Returns
A vector of length getCurSubspaceDim() containing the Ritz values from the previous projected eigensolve.

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

Definition at line 1733 of file AnasaziRTRBase.hpp.

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

Get the index used for extracting Ritz vectors from getRitzVectors().

Because BlockDavidson is a Hermitian solver, all Ritz values are real and all Ritz vectors can be represented in a single column of a multivector. Therefore, getRitzIndex() is not needed when using the output from getRitzVectors().

Returns
An int vector of size getCurSubspaceDim() composed of zeros.

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

Definition at line 1787 of file AnasaziRTRBase.hpp.

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

Get the current residual norms.

Returns
A vector of length getCurSubspaceDim() containing the norms of the residuals, with respect to the orthogonalization manager norm() method.

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

Definition at line 1448 of file AnasaziRTRBase.hpp.

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

Get the current residual 2-norms.

Returns
A vector of length getCurSubspaceDim() containing the 2-norms of the residuals.

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

Definition at line 1462 of file AnasaziRTRBase.hpp.

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

Get the 2-norms of the Ritz residuals.

Returns
A vector of length getCurSubspaceDim() containing the 2-norms of the Ritz residuals.

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

Definition at line 1724 of file AnasaziRTRBase.hpp.

template<class ScalarType , class MV , class OP >
int Anasazi::RTRBase< ScalarType, MV, OP >::getCurSubspaceDim ( ) const
virtual

Get the dimension of the search subspace used to generate the current eigenvectors and eigenvalues.

RTR employs a sequential subspace iteration, maintaining a fixed-rank basis, as opposed to an expanding subspace mechanism employed by Krylov-subspace solvers like BlockKrylovSchur and BlockDavidson.

Returns
An integer specifying the rank of the subspace generated by the eigensolver. If isInitialized() == false, the return is 0. Otherwise, the return will be getBlockSize().

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

Definition at line 1716 of file AnasaziRTRBase.hpp.

template<class ScalarType , class MV , class OP >
int Anasazi::RTRBase< ScalarType, MV, OP >::getMaxSubspaceDim ( ) const
virtual

Get the maximum dimension allocated for the search subspace. For RTR, this always returns getBlockSize().

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

Definition at line 1711 of file AnasaziRTRBase.hpp.

template<class ScalarType , class MV , class OP >
void Anasazi::RTRBase< ScalarType, MV, OP >::setStatusTest ( Teuchos::RCP< StatusTest< ScalarType, MV, OP > >  test)
virtual

Set a new StatusTest for the solver.

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

Definition at line 979 of file AnasaziRTRBase.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > Anasazi::RTRBase< ScalarType, MV, OP >::getStatusTest ( ) const
virtual

Get the current StatusTest used by the solver.

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

Definition at line 989 of file AnasaziRTRBase.hpp.

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

Get a constant reference to the eigenvalue problem.

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

Definition at line 1706 of file AnasaziRTRBase.hpp.

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

Set the blocksize to be used by the iterative solver in solving this eigenproblem.

If the block size is reduced, then the new iterate (and residual and search direction) are chosen as the subset of the current iterate preferred by the sort manager. Otherwise, the solver state is set to uninitialized.

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

Definition at line 704 of file AnasaziRTRBase.hpp.

template<class ScalarType , class MV , class OP >
int Anasazi::RTRBase< ScalarType, MV, OP >::getBlockSize ( ) const
virtual

Get the blocksize to be used by the iterative solver in solving this eigenproblem.

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

Definition at line 1701 of file AnasaziRTRBase.hpp.

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

Set the auxiliary vectors for the solver.

Because the current iterate X cannot be assumed orthogonal to the new auxiliary vectors, a call to setAuxVecs() with a non-empty argument will reset the solver to the uninitialized state.

In order to preserve the current state, the user will need to extract it from the solver using getState(), orthogonalize it against the new auxiliary vectors, and manually reinitialize the solver using initialize().

NOTE: The requirements of the IRTR solvers is such that the auxiliary vectors must be moved into contiguous storage with the current iterate. As a result, the multivector data in auxvecs will be copied, and the multivectors in auxvecs will no longer be referenced. The (unchanged) internal copies of the auxilliary vectors will be made available to the caller by the getAuxVecs() routine. This allows the caller to delete the caller's copies and instead use the copies owned by the solver, avoiding the duplication of data. This is not necessary, however. The partitioning of the auxiliary vectors passed to setAuxVecs() will be preserved.

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

Definition at line 997 of file AnasaziRTRBase.hpp.

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

Get the current auxiliary vectors.

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

Definition at line 1696 of file AnasaziRTRBase.hpp.

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

This method requests that the solver print out its current status to screen.

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

Reimplemented in Anasazi::IRTR< ScalarType, MV, OP >, and Anasazi::SIRTR< ScalarType, MV, OP >.

Definition at line 1593 of file AnasaziRTRBase.hpp.


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