[Trilinos-Users] Anasazi + Eigen values/Eigen Vectors around certain value

Chris Baker cgbaker at gmail.com
Mon Feb 18 08:46:59 MST 2008


Davood,

The standard approach is to use some manner of spectral transformation, so
that the eigenvalues that you want to find are the ones that the solver
wants to find as well.

When looking for interior eigenvalues with a Arnoldi solver (like Anasazi's
BKS), the spectral transformation is typically a shift-invert
transformation. Consider the following generalized eigenvalue problem:
   A*x = B*x*lambda,
where A and B are linear operators, x is a non-zero vector, and lambda is a
scalar.
Remove sigma*B*x (for scalar sigma, your pre-specified value of interest)
from both sides:
  (A - sigma*B)*x = B*x*(lambda - sigma)
Premultiply both sides by the inverse of (A - sigma*B) and 1/(lambda -
sigma):
  inv(A - sigma*B)*B x = x 1/(lambda - sigma)
This is a new (standard) eigenvalue problem, where the eigenvalues closest
to sigma have become the eigenvalues of largest magnitude, i.e., the
eigenvalues which an Arnoldi solve will naturally find:
  Op = inv(A-sigma*B)*B, (the spectral transformation)
  Op*x = x theta,
where
   theta = 1/(lambda-sigma), or
   lambda = 1/theta+sigma   (the back-transformed eigenvalue)
Note also that if A and B are symmetric and B is positive semi-definite,
then the operator Op is symmetric with respect to the inner product defined
by B. In this scenario, you want to pass both Op and B to the
Anasazi::Eigenproblem and specify that the eigenproblem is
symmetric/Hermitian, as in the examples cited below.
Note, in the current version of Anasazi, because Anasazi is not aware of the
spectral transformation present in Op, it cannot back-transform the computed
eigenvalues, i.e., turn the thetas corresponding to the transformed
eigenproblem into the lambdas corresponding to the original eigenproblem.
You must do this after the solver returns. See the examples cited below.

A shift-invert approach requires the ability to "accurately" solve the
shifted system. The
   packages/anasazi/example/BlockKrylovSchur
directory in your Trilinos source tree contains examples for using BKS to
solve shift-invert spectral transformations (with zero shift) via Amesos,
AztecOO or Belos linear solves. Because the shift-invert operator defines a
new eigenvalue problem, the accuracy of the solves determines the
eigenvalues you will get. Poor accuracy in the solvers means that while you
might find accurate thetas for Op*x=x*theta, they will not necessarily
back-transform to accurate lambdas for A*x=B*x*lambda.

It is also possible to use a spectral transformation to use LOBPCG of
BlockDavidson solvers to find interior eigenvalues. These solvers naturally
find the smallest or largest algebraic eigenvalues, i.e., those eigenvalues
furthest to the left or right on the real line (note: they only solve
symmetric/Hermitian eigenvalue problems, so the eigenvalues are all real).
In order to apply one of these solvers, the target eigenvalues must be
transformed to those on the extreme end of the real lines. This can be done
by "folding" the spectrum.

Consider again the eigenvalue problem
  A*x=B*x*lambda
Shift both sides by B*x*sigma:
  (A-sigma*B)*x = B*x*(lambda-sigma)
Premultiply by (A-sigma*B)*inv(B) and you get:
  (A-sigma*B)*inv(B)*(A-sigma*B)*x = B*x*(lambda-sigma)^2
Define the new eigenvalue problem:
  Op = (A-sigma*B)*inv(B)*(A-sigma*B)
  Op*x = B*x*theta,
where
  theta = (lambda-sigma)^2, or
  lambda = sqrt(theta)+sigma (the back-transformed eigenvalue problem).
For this problem, the left-most and right-most eigenvalues are those that
are closest to and furthest away from sigma, respectively. This approach
requires solving against the operator B. If B is the identity (your
eigenproblem is standard), then this solve is trivial. LOBPCG and BD still
may need a preconditioner (for A-sigma*Identity) to perform well.

I don't have much experience with the latter approach, and I can't recommend
any references off the top of my head. Anasazi does not provide any examples
like this.

Rich or Heidi or someone else may have something else to add. A shift-invert
Arnoldi is typically the approach that is taken when finding interior
eigenvalues.

Regards,
Chris

cc: Rich Lehoucq, Heidi Thornquist


On Feb 18, 2008 2:54 AM, Davood Ansari <david.ansari at gmail.com> wrote:

> Hi
>
> I wonder if there is a shortcut approach for finding eigenvalues and
> eigen-vectors
> around a certain pre-specified value (real or complex). I mean using the
> 'which'
> option and just going through the two extreme ends (LM and SM) makes it
> difficult to approach mid spectrum eigen values.
>
> Davood
>
> _______________________________________________
> Trilinos-Users mailing list
> Trilinos-Users at software.sandia.gov
> http://software.sandia.gov/mailman/listinfo/trilinos-users
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: http://software.sandia.gov/mailman/private/trilinos-users/attachments/20080218/252416e9/attachment.html


More information about the Trilinos-Users mailing list