|
Anasazi
Version of the Day
|
This class implements a TraceMIN iteration, a preconditioned iteration for solving linear symmetric positive definite eigenproblems. More...
#include <AnasaziTraceMin.hpp>
Constructor/Destructor | |
| TraceMin (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) | |
| TraceMin constructor with eigenproblem, solver utilities, and parameter list of solver options. | |
This class implements a TraceMIN iteration, a preconditioned iteration for solving linear symmetric positive definite eigenproblems.
TraceMin works by solving a constrained minimization problem
The Y which satisfies this problem is the set of eigenvectors corresponding to the eigenvalues of smallest magnitude. While TraceMin is capable of finding any eigenpairs of
(where K is symmetric and M symmetric positive definite), it converges to the eigenvalues of smallest magnitude of whatever it is given. If a different set of eigenpairs is desired (such as the absolute smallest or the ones closest to a given shift), please perform a spectral transformation before passing the Eigenproblem to the solver.
A TraceMin iteration consists of the following operations:
and
, where
is diagonal.
is an approximation of the eigenpairs of smallest magnitude.The saddle point problem need not be solved to a low relative residual and can be solved either by directly forming the Schur complement or by using a projected Krylov subspace method to solve
Then V is constructed as
. If a preconditioner H is used with the projected Krylov method, it is applied as
and we solve
. (See A Parallel Implementation of the Trace Minimization Eigensolver by Eloy Romero and Jose E. Roman.)
The convergence rate of TraceMin is based on the distribution of eigenvalues. If the eigenvalues are clustered far away from the origin, we have a slow rate of convergence. We can improve our convergence rate by taking advantage of Ritz shifts. Instead of solving
, we solve
.
This method is described in A Trace Minimization Algorithm for the Generalized Eigenvalue Problem, Ahmed H. Sameh, John A. Wisniewski, SIAM Journal on Numerical Analysis, 19(6), pp. 1243-1259 (1982)
Definition at line 114 of file AnasaziTraceMin.hpp.
| Anasazi::Experimental::TraceMin< ScalarType, MV, OP >::TraceMin | ( | 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 | ||
| ) |
TraceMin 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:
"Block Size" - an int specifying the subspace dimension used by the algorithm. This can also be specified using the setBlockSize() method."Maximum Iterations" - an int specifying the maximum number of TraceMin iterations."Saddle Solver Type" - a string specifying the algorithm to use in solving the saddle point problem: "Schur Complement" or "Projected
Krylov". Default: "Projected Krylov""Schur Complement": We explicitly form the Schur complement and use it to solve the linear system. This option may be faster, but it is less numerically stable and does not ensure orthogonality between the current iterate and the update."Projected Krylov": Use a projected iterative method to solve the linear system. If TraceMin was not given a preconditioner, it will use MINRES. Otherwise, it will use GMRES."Shift Type" - a string specifying how to choose Ritz shifts: "No Shift", "Locked Shift", "Trace Leveled Shift", or "Original Shift". Default: "Trace Leveled Shift""No Shift": Do not perform Ritz shifts. This option produces guaranteed convergence but converges linearly. Not recommended."Locked Shift": Do not perform Ritz shifts until an eigenpair is locked. Then, shift based on the largest converged eigenvalue. This option is roughly as safe as "No Shift" but may be faster."Trace Leveled Shift": Do not perform Ritz shifts until the trace of
(i.e. the quantity being minimized has stagnated. Then, shift based on the strategy proposed in The trace minimization method for the symmetric generalized eigenvalue problem. Highly recommended."Original Shift": The original shifting strategy proposed in "The trace minimization method for the symmetric generalized
eigenvalue problem." Compute shifts based on the Ritz values, residuals, and clustering of the eigenvalues. May produce incorrect results for indefinite matrices or small subspace dimensions. Definition at line 196 of file AnasaziTraceMin.hpp.
1.7.6.1