LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ ssygvd()

subroutine ssygvd ( integer  itype,
character  jobz,
character  uplo,
integer  n,
real, dimension( lda, * )  a,
integer  lda,
real, dimension( ldb, * )  b,
integer  ldb,
real, dimension( * )  w,
real, dimension( * )  work,
integer  lwork,
integer, dimension( * )  iwork,
integer  liwork,
integer  info 
)

SSYGVD

Download SSYGVD + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 SSYGVD computes all the eigenvalues, and optionally, the eigenvectors
 of a real generalized symmetric-definite eigenproblem, of the form
 A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.  Here A and
 B are assumed to be symmetric and B is also positive definite.
 If eigenvectors are desired, it uses a divide and conquer algorithm.
Parameters
[in]ITYPE
          ITYPE is INTEGER
          Specifies the problem type to be solved:
          = 1:  A*x = (lambda)*B*x
          = 2:  A*B*x = (lambda)*x
          = 3:  B*A*x = (lambda)*x
[in]JOBZ
          JOBZ is CHARACTER*1
          = 'N':  Compute eigenvalues only;
          = 'V':  Compute eigenvalues and eigenvectors.
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangles of A and B are stored;
          = 'L':  Lower triangles of A and B are stored.
[in]N
          N is INTEGER
          The order of the matrices A and B.  N >= 0.
[in,out]A
          A is REAL array, dimension (LDA, N)
          On entry, the symmetric matrix A.  If UPLO = 'U', the
          leading N-by-N upper triangular part of A contains the
          upper triangular part of the matrix A.  If UPLO = 'L',
          the leading N-by-N lower triangular part of A contains
          the lower triangular part of the matrix A.

          On exit, if JOBZ = 'V', then if INFO = 0, A contains the
          matrix Z of eigenvectors.  The eigenvectors are normalized
          as follows:
          if ITYPE = 1 or 2, Z**T*B*Z = I;
          if ITYPE = 3, Z**T*inv(B)*Z = I.
          If JOBZ = 'N', then on exit the upper triangle (if UPLO='U')
          or the lower triangle (if UPLO='L') of A, including the
          diagonal, is destroyed.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in,out]B
          B is REAL array, dimension (LDB, N)
          On entry, the symmetric matrix B.  If UPLO = 'U', the
          leading N-by-N upper triangular part of B contains the
          upper triangular part of the matrix B.  If UPLO = 'L',
          the leading N-by-N lower triangular part of B contains
          the lower triangular part of the matrix B.

          On exit, if INFO <= N, the part of B containing the matrix is
          overwritten by the triangular factor U or L from the Cholesky
          factorization B = U**T*U or B = L*L**T.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[out]W
          W is REAL array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]WORK
          WORK is REAL array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.
          If N <= 1,               LWORK >= 1.
          If JOBZ = 'N' and N > 1, LWORK >= 2*N+1.
          If JOBZ = 'V' and N > 1, LWORK >= 1 + 6*N + 2*N**2.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal sizes of the WORK and IWORK
          arrays, returns these values as the first entries of the WORK
          and IWORK arrays, and no error message related to LWORK or
          LIWORK is issued by XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of the array IWORK.
          If N <= 1,                LIWORK >= 1.
          If JOBZ  = 'N' and N > 1, LIWORK >= 1.
          If JOBZ  = 'V' and N > 1, LIWORK >= 3 + 5*N.

          If LIWORK = -1, then a workspace query is assumed; the
          routine only calculates the optimal sizes of the WORK and
          IWORK arrays, returns these values as the first entries of
          the WORK and IWORK arrays, and no error message related to
          LWORK or LIWORK is issued by XERBLA.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  SPOTRF or SSYEVD returned an error code:
             <= N:  if INFO = i and JOBZ = 'N', then the algorithm
                    failed to converge; i off-diagonal elements of an
                    intermediate tridiagonal form did not converge to
                    zero;
                    if INFO = i and JOBZ = 'V', then the algorithm
                    failed to compute an eigenvalue while working on
                    the submatrix lying in rows and columns INFO/(N+1)
                    through mod(INFO,N+1);
             > N:   if INFO = N + i, for 1 <= i <= N, then the leading
                    principal minor of order i of B is not positive.
                    The factorization of B could not be completed and
                    no eigenvalues or eigenvectors were computed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
  Modified so that no backsubstitution is performed if SSYEVD fails to
  converge (NEIG in old code could be greater than N causing out of
  bounds reference to A - reported by Ralf Meyer).  Also corrected the
  description of INFO and the test on ITYPE. Sven, 16 Feb 05.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 219 of file ssygvd.f.

221*
222* -- LAPACK driver routine --
223* -- LAPACK is a software package provided by Univ. of Tennessee, --
224* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
225*
226* .. Scalar Arguments ..
227 CHARACTER JOBZ, UPLO
228 INTEGER INFO, ITYPE, LDA, LDB, LIWORK, LWORK, N
229* ..
230* .. Array Arguments ..
231 INTEGER IWORK( * )
232 REAL A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
233* ..
234*
235* =====================================================================
236*
237* .. Parameters ..
238 REAL ONE
239 parameter( one = 1.0e+0 )
240* ..
241* .. Local Scalars ..
242 LOGICAL LQUERY, UPPER, WANTZ
243 CHARACTER TRANS
244 INTEGER LIOPT, LIWMIN, LOPT, LWMIN
245* ..
246* .. External Functions ..
247 LOGICAL LSAME
248 REAL SROUNDUP_LWORK
249 EXTERNAL lsame, sroundup_lwork
250* ..
251* .. External Subroutines ..
252 EXTERNAL spotrf, ssyevd, ssygst, strmm, strsm, xerbla
253* ..
254* .. Intrinsic Functions ..
255 INTRINSIC max, real
256* ..
257* .. Executable Statements ..
258*
259* Test the input parameters.
260*
261 wantz = lsame( jobz, 'V' )
262 upper = lsame( uplo, 'U' )
263 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
264*
265 info = 0
266 IF( n.LE.1 ) THEN
267 liwmin = 1
268 lwmin = 1
269 ELSE IF( wantz ) THEN
270 liwmin = 3 + 5*n
271 lwmin = 1 + 6*n + 2*n**2
272 ELSE
273 liwmin = 1
274 lwmin = 2*n + 1
275 END IF
276 lopt = lwmin
277 liopt = liwmin
278 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
279 info = -1
280 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
281 info = -2
282 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
283 info = -3
284 ELSE IF( n.LT.0 ) THEN
285 info = -4
286 ELSE IF( lda.LT.max( 1, n ) ) THEN
287 info = -6
288 ELSE IF( ldb.LT.max( 1, n ) ) THEN
289 info = -8
290 END IF
291*
292 IF( info.EQ.0 ) THEN
293 work( 1 ) = sroundup_lwork(lopt)
294 iwork( 1 ) = liopt
295*
296 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
297 info = -11
298 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
299 info = -13
300 END IF
301 END IF
302*
303 IF( info.NE.0 ) THEN
304 CALL xerbla( 'SSYGVD', -info )
305 RETURN
306 ELSE IF( lquery ) THEN
307 RETURN
308 END IF
309*
310* Quick return if possible
311*
312 IF( n.EQ.0 )
313 $ RETURN
314*
315* Form a Cholesky factorization of B.
316*
317 CALL spotrf( uplo, n, b, ldb, info )
318 IF( info.NE.0 ) THEN
319 info = n + info
320 RETURN
321 END IF
322*
323* Transform problem to standard eigenvalue problem and solve.
324*
325 CALL ssygst( itype, uplo, n, a, lda, b, ldb, info )
326 CALL ssyevd( jobz, uplo, n, a, lda, w, work, lwork, iwork, liwork,
327 $ info )
328 lopt = int( max( real( lopt ), real( work( 1 ) ) ) )
329 liopt = int( max( real( liopt ), real( iwork( 1 ) ) ) )
330*
331 IF( wantz .AND. info.EQ.0 ) THEN
332*
333* Backtransform eigenvectors to the original problem.
334*
335 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
336*
337* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
338* backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
339*
340 IF( upper ) THEN
341 trans = 'N'
342 ELSE
343 trans = 'T'
344 END IF
345*
346 CALL strsm( 'Left', uplo, trans, 'Non-unit', n, n, one,
347 $ b, ldb, a, lda )
348*
349 ELSE IF( itype.EQ.3 ) THEN
350*
351* For B*A*x=(lambda)*x;
352* backtransform eigenvectors: x = L*y or U**T*y
353*
354 IF( upper ) THEN
355 trans = 'T'
356 ELSE
357 trans = 'N'
358 END IF
359*
360 CALL strmm( 'Left', uplo, trans, 'Non-unit', n, n, one,
361 $ b, ldb, a, lda )
362 END IF
363 END IF
364*
365 work( 1 ) = sroundup_lwork(lopt)
366 iwork( 1 ) = liopt
367*
368 RETURN
369*
370* End of SSYGVD
371*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ssyevd(jobz, uplo, n, a, lda, w, work, lwork, iwork, liwork, info)
SSYEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices
Definition ssyevd.f:176
subroutine ssygst(itype, uplo, n, a, lda, b, ldb, info)
SSYGST
Definition ssygst.f:127
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine spotrf(uplo, n, a, lda, info)
SPOTRF
Definition spotrf.f:107
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine strmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
STRMM
Definition strmm.f:177
subroutine strsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
STRSM
Definition strsm.f:181
Here is the call graph for this function:
Here is the caller graph for this function: