Scilab Function
Last update : 25/10/2006

dsaupd - Interface for the Implicitly Restarted Arnoldi Iteration, to compute approximations to a few eigenpairs of a real and symmetric linear operator

Calling Sequence

[IDO,RESID,V,IPARAM,IPNTR,WORKD,WORKL,INFO] = dsaupd(ID0,BMAT,N,WHICH,NEV,TOL,RESID,NCV,V,IPARAM,IPNTR,WORKD,WORKL,INFO)

Parameters

Description

Reverse communication interface for the Implicitly Restarted Arnoldi Iteration. For symmetric problems this reduces to a variant of the Lanczos method. This method has been designed to compute approximations to a few eigenpairs of a linear operator OP that is real and symmetric with respect to a real positive semi-definite symmetric matrix B, i.e. B*OP = (OP`)*B .

Another way to express this condition is < x,OPy > = < OPx,y > where <z,w > = z`Bw .

In the standard eigenproblem B is the identity matrix.( A` denotes transpose of A)

The computed approximate eigenvalues are called Ritz values and the corresponding approximate eigenvectors are called Ritz vectors.

dsaupd is usually called iteratively to solve one of the following problems:

Mode 1: A*x = lambda*x, A symmetric ===> OP = A and B = I.

Mode 2: A*x = lambda*M*x, A symmetric, M symmetric positive definite ===> OP = inv[M]*A and B = M. ===> (If M can be factored see remark 3 below)

Mode 3: K*x = lambda*M*x, K symmetric, M symmetric positive semi-definite ===> OP = (inv[K - sigma*M])*M and B = M. ===> Shift-and-Invert mode

Mode 4: K*x = lambda*KG*x, K symmetric positive semi-definite, KG symmetric indefinite ===> OP = (inv[K - sigma*KG])*K and B = K. ===> Buckling mode

Mode 5: A*x = lambda*M*x, A symmetric, M symmetric positive semi-definite ===> OP = inv[A - sigma*M]*[A + sigma*M] and B = M. ===> Cayley transformed mode

NOTE: The action of w <- inv[A - sigma*M]*v or w <- inv[M]*v should be accomplished either by a direct method using a sparse matrix factorization and solving [A - sigma*M]*w = v or M*w = v ,

or through an iterative method for solving these systems. If an iterative method is used, the convergence test must be more stringent than the accuracy requirements for the eigenvalue approximations.

Remarks

1. The converged Ritz values are always returned in ascending algebraic order. The computed Ritz values are approximate eigenvalues of OP. The selection of WHICH should be made with this in mind when Mode = 3,4,5. After convergence, approximate eigenvalues of the original problem may be obtained with the ARPACK subroutine dseupd .

2. If the Ritz vectors corresponding to the converged Ritz values are needed, the user must call dseupd immediately following completion of dsaupd . This is new starting with version 2.1 of ARPACK.

3. If M can be factored into a Cholesky factorization M = LL` then Mode = 2 should not be selected. Instead one should use Mode = 1 with OP = inv(L)*A*inv(L`). Appropriate triangular linear systems should be solved with L and L` rather than computing inverses. After convergence, an approximate eigenvector z of the original problem is recovered by solving L`z = x where x is a Ritz vector of OP.

4. At present there is no a-priori analysis to guide the selection of NCV relative to NEV. The only formal requrement is that NCV > NEV. However, it is recommended that NCV >= 2*NEV. If many problems of the same type are to be solved, one should experiment with increasing NCV while keeping NEV fixed for a given test problem. This will usually decrease the required number of OP*x operations but it also increases the work and storage required to maintain the orthogonal basis vectors. The optimal "cross-over" with respect to CPU time is problem dependent and must be determined empirically.

5. If IPARAM(7) = 2 then in the Reverse commuication interface the user must do the following. When IDO = 1, Y = OP * X is to be computed. When IPARAM(7) = 2 OP = inv(B)*A. After computing A*X the user must overwrite X with A*X. Y is then the solution to the linear set of equations B*Y = A*X.

6. When IPARAM(1) = 0, and IDO = 3, the user needs to provide the NP = IPARAM(8) shifts in locations: 1 WORKL(IPNTR(11)) 2 WORKL(IPNTR(11)+1) NP WORKL(IPNTR(11)+NP-1). The eigenvalues of the current tridiagonal matrix are located in WORKL(IPNTR(6)) through WORKL(IPNTR(6)+NCV-1). They are in the order defined by WHICH. The associated Ritz estimates are located in WORKL(IPNTR(8)), WORKL(IPNTR(8)+1), ... , WORKL(IPNTR(8)+NCV-1).

See Also

dnaupd ,  

Authors

Danny Sorensen, Richard Lehoucq, Phuong Vu CRPC / Rice University Applied Mathematics Rice University Houston, Texas

Bibliography

1. D.C. Sorensen, "Implicit Application of Polynomial Filters in a k-Step Arnoldi Method", SIAM J. Matr. Anal. Apps., 13 (1992), pp 357-385.

2. R.B. Lehoucq, "Analysis and Implementation of an Implicitly Restarted Arnoldi Iteration", Rice University Technical Report TR95-13, Department of Computational and Applied Mathematics.

3. B.N. Parlett and Y. Saad, "Complex Shift and Invert Strategies for Real Matrices", Linear Algebra and its Applications, vol 88/89, pp 575-595, (1987).

Used Function

Based on ARPACK routine dnaupd