Anasazi
Version of the Day
|
This is an abstract base class for the trace minimization eigensolvers. More...
#include <AnasaziTraceMinBase.hpp>
Public Member Functions | |
Constructor/Destructor | |
TraceMinBase (const RCP< Eigenproblem< ScalarType, MV, OP > > &problem, const RCP< SortManager< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > > &sorter, const RCP< OutputManager< ScalarType > > &printer, const RCP< StatusTest< ScalarType, MV, OP > > &tester, const RCP< MatOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList ¶ms) | |
TraceMinBase constructor with eigenproblem, solver utilities, and parameter list of solver options. More... | |
virtual | ~TraceMinBase () |
Anasazi::TraceMinBase destructor. More... | |
Solver methods | |
void | iterate () |
This method performs trace minimization iterations until the status test indicates the need to stop or an error occurs (in which case, an appropriate exception is thrown). More... | |
void | harmonicIterate () |
void | initialize (TraceMinBaseState< ScalarType, MV > &newstate) |
Initialize the solver to an iterate, optionally providing the other members of the state. More... | |
void | harmonicInitialize (TraceMinBaseState< ScalarType, MV > newstate) |
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... | |
TraceMinBaseState< 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... | |
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 the trace minimization methods, this always returns numBlocks*blockSize. More... | |
Accessor routines from Eigensolver | |
void | setStatusTest (RCP< StatusTest< ScalarType, MV, OP > > test) |
Set a new StatusTest for the solver. More... | |
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< RCP< const MV > > &auxvecs) |
Set the auxiliary vectors for the solver. More... | |
Teuchos::Array< RCP< const MV > > | getAuxVecs () const |
Get the auxiliary vectors for the solver. More... | |
BlockBase-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 is an abstract base class for the trace minimization eigensolvers.
For more information, please see Anasazi::TraceMin (with constant subspace dimension) and Anasazi::TraceMinDavidson (with expanding subspaces)
Definition at line 150 of file AnasaziTraceMinBase.hpp.
Anasazi::Experimental::TraceMinBase< ScalarType, MV, OP >::TraceMinBase | ( | const RCP< Eigenproblem< ScalarType, MV, OP > > & | problem, |
const RCP< SortManager< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > > & | sorter, | ||
const RCP< OutputManager< ScalarType > > & | printer, | ||
const RCP< StatusTest< ScalarType, MV, OP > > & | tester, | ||
const RCP< MatOrthoManager< ScalarType, MV, OP > > & | ortho, | ||
Teuchos::ParameterList & | params | ||
) |
TraceMinBase 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:
"Saddle Solver Type"
- a string
specifying how to solve the saddle point problem arising at each iteration. Options are "Projected Krylov", "Schur Complement", and "Block Diagonal Preconditioned Minres". Default: "Projected Krylov""Projected Krylov"
: Uses projected-minres to solve the problem."Schur Complement"
: Explicitly forms the (inexact) Schur complement using minres."Block Diagonal Preconditioned Minres"
: Uses a block preconditioner on the entire saddle point problem. For more information, please see "Overview of Anasazi and its newest eigensolver, TraceMin" on the main Anasazi page. We recommend using "Projected Krylov" in the absence of preconditioning. If you want to use a preconditioner, "Block Diagonal Preconditioned Minres" is recommended. "Schur Complement" mainly exists for special use cases."When To Shift"
- a string
specifying when Ritz shifts should be performed. Options are "Never", "After Trace Levels", and "Always". Default: "Always""Never"
: Do not perform Ritz shifts. This option produces guaranteed convergence but converges linearly. Not recommended."After Trace Levels"
: Do not perform Ritz shifts until the trace of has stagnated (i.e. the relative change in trace has become small). The MagnitudeType
specifying how small the relative change in trace must become may be provided via the parameter "Trace Threshold"
, whose default value is 0.02."Always"
: Always attempt to use Ritz shifts."How To Choose Shift"
- a string
specifying how to choose the Ritz shifts (assuming Ritz shifts are being used). Options are "Largest Converged", "Adjusted Ritz Values", and "Ritz Values". Default: "Adjusted Ritz Values""Largest Converged"
: Ritz shifts are chosen to be the largest converged eigenvalue. Until an eigenvalue converges, the Ritz shifts are all 0."Adjusted Ritz Values"
: Ritz shifts are chosen based on the Ritz values and their associated residuals in such a way as to guarantee global convergence. This method is described in "The trace minimization method for the symmetric generalized eigenvalue problem.""Ritz Values"
: Ritz shifts are chosen to equal the Ritz values. This does NOT guarantee global convergence."Use Multiple Shifts"
- a bool
specifying whether to use one or many Ritz shifts (assuming shifting is enabled). Default: trueAnasazi's trace minimization solvers are still in development, and we plan to add additional features in the future including additional saddle point solvers.
Definition at line 621 of file AnasaziTraceMinBase.hpp.
|
virtual |
Anasazi::TraceMinBase destructor.
Definition at line 745 of file AnasaziTraceMinBase.hpp.
|
virtual |
This method performs trace minimization 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 initialized; if not, it will call initialize(). After initialization, the solver performs TraceMin iterations until the status test evaluates as Passed, at which point the method returns to the caller.
The trace minimization 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 TraceMinBase-specific exceptions.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 897 of file AnasaziTraceMinBase.hpp.
void Anasazi::Experimental::TraceMinBase< ScalarType, MV, OP >::initialize | ( | TraceMinBaseState< ScalarType, MV > & | newstate | ) |
Initialize the solver to an iterate, optionally providing the other members of the state.
The TraceMinBase 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.
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 1135 of file AnasaziTraceMinBase.hpp.
|
virtual |
Initialize the solver with the initial vectors from the eigenproblem or random data.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 1114 of file AnasaziTraceMinBase.hpp.
|
virtual |
Indicates whether the solver has been initialized or not.
true:
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 2269 of file AnasaziTraceMinBase.hpp.
TraceMinBaseState< ScalarType, MV > Anasazi::Experimental::TraceMinBase< ScalarType, MV, OP >::getState | ( | ) | const |
Get access to the current state of the eigensolver.
The data is only valid if isInitialized() == true
.
Definition at line 857 of file AnasaziTraceMinBase.hpp.
|
virtual |
Get the current iteration count.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 849 of file AnasaziTraceMinBase.hpp.
|
virtual |
Reset the iteration count.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 841 of file AnasaziTraceMinBase.hpp.
|
virtual |
Get access to the current Ritz vectors.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 833 of file AnasaziTraceMinBase.hpp.
|
virtual |
Get the Ritz values for the previous iteration.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 820 of file AnasaziTraceMinBase.hpp.
|
virtual |
Get the index used for extracting individual Ritz vectors from getRitzVectors().
Because the trace minimization methods are a Hermitian solvers, 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 811 of file AnasaziTraceMinBase.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 2465 of file AnasaziTraceMinBase.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 2507 of file AnasaziTraceMinBase.hpp.
|
virtual |
Get the 2-norms of the residuals.
The Ritz residuals are not defined for trace minimization iterations. 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 803 of file AnasaziTraceMinBase.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 793 of file AnasaziTraceMinBase.hpp.
|
virtual |
Get the maximum dimension allocated for the search subspace. For the trace minimization methods, this always returns numBlocks*blockSize.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 786 of file AnasaziTraceMinBase.hpp.
|
virtual |
Set a new StatusTest for the solver.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 2549 of file AnasaziTraceMinBase.hpp.
|
virtual |
Get the current StatusTest used by the solver.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 2558 of file AnasaziTraceMinBase.hpp.
|
virtual |
Get a constant reference to the eigenvalue problem.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 778 of file AnasaziTraceMinBase.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 TraceMinBase 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 subspace dimension at each iteration.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 752 of file AnasaziTraceMinBase.hpp.
|
virtual |
Get the blocksize used by the iterative solver.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 770 of file AnasaziTraceMinBase.hpp.
|
virtual |
Set the auxiliary vectors for the solver.
Auxiliary vectors are ones that you want your eigenvectors to be held orthogonal to. One example of where you may want to use this is in the computation of the Fiedler vector, where you would likely want to project against the vector of all 1s.
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 2421 of file AnasaziTraceMinBase.hpp.
|
virtual |
Get the auxiliary vectors for the solver.
Implements Anasazi::Eigensolver< ScalarType, MV, OP >.
Definition at line 762 of file AnasaziTraceMinBase.hpp.
void Anasazi::Experimental::TraceMinBase< 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.
Definition at line 2275 of file AnasaziTraceMinBase.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 2566 of file AnasaziTraceMinBase.hpp.