Anasazi
Version of the Day
|
This class implements a Block Davidson iteration, a preconditioned iteration for solving linear Hermitian eigenproblems. More...
#include <AnasaziBlockDavidson.hpp>
Public Member Functions | |
Constructor/Destructor | |
BlockDavidson (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< MatOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList ¶ms) | |
BlockDavidson constructor with eigenproblem, solver utilities, and parameter list of solver options. More... | |
virtual | ~BlockDavidson () |
Anasazi::BlockDavidson destructor. More... | |
Solver methods | |
void | iterate () |
This method performs BlockDavidson iterations until the status test indicates the need to stop or an error occurs (in which case, an appropriate exception is thrown). More... | |
void | initialize (BlockDavidsonState< ScalarType, MV > &newstate) |
Initialize the solver to an iterate, optionally providing the current basis and projected problem matrix, the current Ritz vectors and values, and the current 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... | |
BlockDavidsonState< ScalarType, MV > | getState () const |
Get access to 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 access to the current Ritz vectors. More... | |
std::vector< Value< ScalarType > > | getRitzValues () |
Get the Ritz values for the previous iteration. More... | |
std::vector< int > | getRitzIndex () |
Get the index used for extracting individual Ritz vectors from getRitzVectors(). More... | |
std::vector< typename Teuchos::ScalarTraits < ScalarType >::magnitudeType > | getResNorms () |
Get the current residual norms, computing the norms if they are not up-to-date with the current residual vectors. More... | |
std::vector< typename Teuchos::ScalarTraits < ScalarType >::magnitudeType > | getRes2Norms () |
Get the current residual 2-norms, computing the norms if they are not up-to-date with the current residual vectors. More... | |
std::vector< typename Teuchos::ScalarTraits < ScalarType >::magnitudeType > | getRitzRes2Norms () |
Get the 2-norms of the 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 BlockDavidson, this always returns numBlocks*blockSize. 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. More... | |
int | getBlockSize () const |
Get the blocksize used by the iterative solver. 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 auxiliary vectors for the solver. More... | |
BlockDavidson-specific accessor routines | |
void | setSize (int blockSize, int numBlocks) |
Set the blocksize and number of blocks to be used by the iterative solver in solving this eigenproblem. More... | |
Output methods | |
void | currentStatus (std::ostream &os) |
This method requests that the solver print out its current status to the given output stream. 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 ¶ms) | |
Basic Constructor. More... | |
virtual | ~Eigensolver () |
Destructor. More... | |
This class implements a Block Davidson iteration, a preconditioned iteration for solving linear Hermitian eigenproblems.
This method is described in A Comparison of Eigensolvers for Large-scale 3D Modal Analysis Using AMG-Preconditioned Iterative Methods, P. Arbenz, U. L. Hetmaniuk, R. B. Lehoucq, R. S. Tuminaro, Internat. J. for Numer. Methods Engrg., 64, pp. 204-236 (2005)
Definition at line 164 of file AnasaziBlockDavidson.hpp.
Anasazi::BlockDavidson< ScalarType, MV, OP >::BlockDavidson | ( | 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< MatOrthoManager< ScalarType, MV, OP > > & | ortho, | ||
Teuchos::ParameterList & | params | ||
) |
BlockDavidson constructor with eigenproblem, solver utilities, and parameter list of solver options.
This constructor takes pointers required by the eigensolver, in addition to a parameter list of options for the eigensolver. These options include the following:
int
specifying the block size used by the algorithm. This can also be specified using the setBlockSize() method.int
specifying the maximum number of blocks allocated for the solver basis. Definition at line 553 of file AnasaziBlockDavidson.hpp.
|
virtual |
Anasazi::BlockDavidson destructor.
Definition at line 629 of file AnasaziBlockDavidson.hpp.
|
virtual |
This method performs BlockDavidson iterations until the status test indicates the need to stop or an error occurs (in which case, an appropriate exception is thrown).
iterate() will first determine whether the solver is uninitialized; if not, it will call initialize(). After initialization, the solver performs block Davidson iterations until the status test evaluates as Passed, at which point the method returns to the caller.
The block Davidson iteration proceeds as follows:
The status test is queried at the beginning of the iteration.
Possible exceptions thrown include std::invalid_argument or one of the BlockDavidson-specific exceptions.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 1277 of file AnasaziBlockDavidson.hpp.
void Anasazi::BlockDavidson< ScalarType, MV, OP >::initialize | ( | BlockDavidsonState< ScalarType, MV > & | newstate | ) |
Initialize the solver to an iterate, optionally providing the current basis and projected problem matrix, the current Ritz vectors and values, and the current residual.
The BlockDavidson eigensolver contains a certain amount of state, including the current Krylov basis, the current eigenvectors, the current residual, etc. (see getState())
initialize() gives the user the opportunity to manually set these, although this must be done with caution, as the validity of the user input will not be checked.
Only the first newstate.curDim
columns of newstate.V
and newstate.KK
and the first newstate.curDim
rows of newstate.KK
will be used.
If newstate.V == getState().V
, then the data is not copied. The same holds for newstate.KK
, newstate.X
, newstate.KX
, newstate.MX
, and newstate.R
Only the upper triangular half of newstate.KK
is used to initialize the state of the solver.
true
(see post-conditions of isInitialize())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., KX) not given to initialize() will be generated.
Note, for any pointer in newstate
which directly points to the multivectors in the solver, the data is not copied.
Definition at line 898 of file AnasaziBlockDavidson.hpp.
|
virtual |
Initialize the solver with the initial vectors from the eigenproblem or random data.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 1267 of file AnasaziBlockDavidson.hpp.
|
virtual |
Indicates whether the solver has been initialized or not.
true:
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 768 of file AnasaziBlockDavidson.hpp.
BlockDavidsonState< ScalarType, MV > Anasazi::BlockDavidson< ScalarType, MV, OP >::getState | ( | ) | const |
Get access to the current state of the eigensolver.
The data is only valid if isInitialized() == true
.
The data for the preconditioned residual is only meaningful in the scenario that the solver throws a ::BlockDavidsonRitzFailure exception during iterate().
Definition at line 741 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get the current iteration count.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 733 of file AnasaziBlockDavidson.hpp.
|
virtual |
Reset the iteration count.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 725 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get access to the current Ritz vectors.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 717 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get the Ritz values for the previous iteration.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 704 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get the index used for extracting individual 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().
int
vector of size getCurSubspaceDim() composed of zeros. Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 695 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get the current residual norms, computing the norms if they are not up-to-date with the current residual vectors.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 1536 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get the current residual 2-norms, computing the norms if they are not up-to-date with the current residual vectors.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 1550 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get the 2-norms of the residuals.
The Ritz residuals are not defined for the LOBPCG iteration. Hence, this method returns the 2-norms of the direct residuals, and is equivalent to calling getRes2Norms().
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 687 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get the dimension of the search subspace used to generate the current eigenvectors and eigenvalues.
false
, the return is 0. Otherwise, it will be some strictly positive multiple of getBlockSize(). Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 677 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get the maximum dimension allocated for the search subspace. For BlockDavidson, this always returns numBlocks*blockSize.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 670 of file AnasaziBlockDavidson.hpp.
|
virtual |
Set a new StatusTest for the solver.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 1562 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get the current StatusTest used by the solver.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 1571 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get a constant reference to the eigenvalue problem.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 662 of file AnasaziBlockDavidson.hpp.
|
virtual |
Set the blocksize.
This method is required to support the interface provided by Eigensolver. However, the preferred method of setting the allocated size for the BlockDavidson eigensolver is setSize(). In fact, setBlockSize() simply calls setSize(), maintaining the current number of blocks.
The block size determines the number of Ritz vectors and values that are computed on each iteration, thereby determining the increase in the Krylov subspace at each iteration.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 636 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get the blocksize used by the iterative solver.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 654 of file AnasaziBlockDavidson.hpp.
|
virtual |
Set the auxiliary vectors for the solver.
Because the current basis V cannot be assumed orthogonal to the new auxiliary vectors, a call to setAuxVecs() will reset the solver to the uninitialized state. This happens only in the case where the new auxiliary vectors have a combined dimension of greater than zero.
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 reinitialize using initialize().
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 861 of file AnasaziBlockDavidson.hpp.
|
virtual |
Get the auxiliary vectors for the solver.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 645 of file AnasaziBlockDavidson.hpp.
void Anasazi::BlockDavidson< ScalarType, MV, OP >::setSize | ( | int | blockSize, |
int | numBlocks | ||
) |
Set the blocksize and number of blocks to be used by the iterative solver in solving this eigenproblem.
Changing either the block size or the number of blocks will reset the solver to an uninitialized state.
The requested block size must be strictly positive; the number of blocks must be greater than one. Invalid arguments will result in a std::invalid_argument exception.
Definition at line 774 of file AnasaziBlockDavidson.hpp.
|
virtual |
This method requests that the solver print out its current status to the given output stream.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 1723 of file AnasaziBlockDavidson.hpp.