LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine chegvd ( integer  ITYPE,
character  JOBZ,
character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldb, * )  B,
integer  LDB,
real, dimension( * )  W,
complex, dimension( * )  WORK,
integer  LWORK,
real, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

CHEGVD

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

Purpose:
 CHEGVD computes all the eigenvalues, and optionally, the eigenvectors
 of a complex generalized Hermitian-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 Hermitian and B is also positive definite.
 If eigenvectors are desired, it uses a divide and conquer algorithm.

 The divide and conquer algorithm makes very mild assumptions about
 floating point arithmetic. It will work on machines with a guard
 digit in add/subtract, or on those binary machines without guard
 digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
 Cray-2. It could conceivably fail on hexadecimal or decimal machines
 without guard digits, but we know of none.
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 COMPLEX array, dimension (LDA, N)
          On entry, the Hermitian 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**H*B*Z = I;
          if ITYPE = 3, Z**H*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 COMPLEX array, dimension (LDB, N)
          On entry, the Hermitian 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**H*U or B = L*L**H.
[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 COMPLEX array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The length of the array WORK.
          If N <= 1,                LWORK >= 1.
          If JOBZ  = 'N' and N > 1, LWORK >= N + 1.
          If JOBZ  = 'V' and N > 1, LWORK >= 2*N + N**2.

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

          If LRWORK = -1, then a workspace query is assumed; the
          routine only calculates the optimal sizes of the WORK, RWORK
          and IWORK arrays, returns these values as the first entries
          of the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK 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, RWORK
          and IWORK arrays, returns these values as the first entries
          of the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK 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:  CPOTRF or CHEEVD 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
                    minor of order i of B is not positive definite.
                    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.
Date
November 2015
Further Details:
  Modified so that no backsubstitution is performed if CHEEVD 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 251 of file chegvd.f.

251 *
252 * -- LAPACK driver routine (version 3.6.0) --
253 * -- LAPACK is a software package provided by Univ. of Tennessee, --
254 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
255 * November 2015
256 *
257 * .. Scalar Arguments ..
258  CHARACTER jobz, uplo
259  INTEGER info, itype, lda, ldb, liwork, lrwork, lwork, n
260 * ..
261 * .. Array Arguments ..
262  INTEGER iwork( * )
263  REAL rwork( * ), w( * )
264  COMPLEX a( lda, * ), b( ldb, * ), work( * )
265 * ..
266 *
267 * =====================================================================
268 *
269 * .. Parameters ..
270  COMPLEX cone
271  parameter ( cone = ( 1.0e+0, 0.0e+0 ) )
272 * ..
273 * .. Local Scalars ..
274  LOGICAL lquery, upper, wantz
275  CHARACTER trans
276  INTEGER liopt, liwmin, lopt, lropt, lrwmin, lwmin
277 * ..
278 * .. External Functions ..
279  LOGICAL lsame
280  EXTERNAL lsame
281 * ..
282 * .. External Subroutines ..
283  EXTERNAL cheevd, chegst, cpotrf, ctrmm, ctrsm, xerbla
284 * ..
285 * .. Intrinsic Functions ..
286  INTRINSIC max, real
287 * ..
288 * .. Executable Statements ..
289 *
290 * Test the input parameters.
291 *
292  wantz = lsame( jobz, 'V' )
293  upper = lsame( uplo, 'U' )
294  lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
295 *
296  info = 0
297  IF( n.LE.1 ) THEN
298  lwmin = 1
299  lrwmin = 1
300  liwmin = 1
301  ELSE IF( wantz ) THEN
302  lwmin = 2*n + n*n
303  lrwmin = 1 + 5*n + 2*n*n
304  liwmin = 3 + 5*n
305  ELSE
306  lwmin = n + 1
307  lrwmin = n
308  liwmin = 1
309  END IF
310  lopt = lwmin
311  lropt = lrwmin
312  liopt = liwmin
313  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
314  info = -1
315  ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
316  info = -2
317  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
318  info = -3
319  ELSE IF( n.LT.0 ) THEN
320  info = -4
321  ELSE IF( lda.LT.max( 1, n ) ) THEN
322  info = -6
323  ELSE IF( ldb.LT.max( 1, n ) ) THEN
324  info = -8
325  END IF
326 *
327  IF( info.EQ.0 ) THEN
328  work( 1 ) = lopt
329  rwork( 1 ) = lropt
330  iwork( 1 ) = liopt
331 *
332  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
333  info = -11
334  ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
335  info = -13
336  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
337  info = -15
338  END IF
339  END IF
340 *
341  IF( info.NE.0 ) THEN
342  CALL xerbla( 'CHEGVD', -info )
343  RETURN
344  ELSE IF( lquery ) THEN
345  RETURN
346  END IF
347 *
348 * Quick return if possible
349 *
350  IF( n.EQ.0 )
351  $ RETURN
352 *
353 * Form a Cholesky factorization of B.
354 *
355  CALL cpotrf( uplo, n, b, ldb, info )
356  IF( info.NE.0 ) THEN
357  info = n + info
358  RETURN
359  END IF
360 *
361 * Transform problem to standard eigenvalue problem and solve.
362 *
363  CALL chegst( itype, uplo, n, a, lda, b, ldb, info )
364  CALL cheevd( jobz, uplo, n, a, lda, w, work, lwork, rwork, lrwork,
365  $ iwork, liwork, info )
366  lopt = max( REAL( LOPT ), REAL( WORK( 1 ) ) )
367  lropt = max( REAL( LROPT ), REAL( RWORK( 1 ) ) )
368  liopt = max( REAL( LIOPT ), REAL( IWORK( 1 ) ) )
369 *
370  IF( wantz .AND. info.EQ.0 ) THEN
371 *
372 * Backtransform eigenvectors to the original problem.
373 *
374  IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
375 *
376 * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
377 * backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
378 *
379  IF( upper ) THEN
380  trans = 'N'
381  ELSE
382  trans = 'C'
383  END IF
384 *
385  CALL ctrsm( 'Left', uplo, trans, 'Non-unit', n, n, cone,
386  $ b, ldb, a, lda )
387 *
388  ELSE IF( itype.EQ.3 ) THEN
389 *
390 * For B*A*x=(lambda)*x;
391 * backtransform eigenvectors: x = L*y or U**H *y
392 *
393  IF( upper ) THEN
394  trans = 'C'
395  ELSE
396  trans = 'N'
397  END IF
398 *
399  CALL ctrmm( 'Left', uplo, trans, 'Non-unit', n, n, cone,
400  $ b, ldb, a, lda )
401  END IF
402  END IF
403 *
404  work( 1 ) = lopt
405  rwork( 1 ) = lropt
406  iwork( 1 ) = liopt
407 *
408  RETURN
409 *
410 * End of CHEGVD
411 *
subroutine chegst(ITYPE, UPLO, N, A, LDA, B, LDB, INFO)
CHEGST
Definition: chegst.f:129
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ctrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRSM
Definition: ctrsm.f:182
subroutine cheevd(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
CHEEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices ...
Definition: cheevd.f:207
subroutine cpotrf(UPLO, N, A, LDA, INFO)
CPOTRF
Definition: cpotrf.f:109
subroutine ctrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
CTRMM
Definition: ctrmm.f:179
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: