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

◆ zhegvd()

subroutine zhegvd ( integer  itype,
character  jobz,
character  uplo,
integer  n,
complex*16, dimension( lda, * )  a,
integer  lda,
complex*16, dimension( ldb, * )  b,
integer  ldb,
double precision, dimension( * )  w,
complex*16, dimension( * )  work,
integer  lwork,
double precision, dimension( * )  rwork,
integer  lrwork,
integer, dimension( * )  iwork,
integer  liwork,
integer  info 
)

ZHEGVD

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

Purpose:
 ZHEGVD 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.
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*16 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*16 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 DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]WORK
          WORK is COMPLEX*16 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 DOUBLE PRECISION 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:  ZPOTRF or ZHEEVD 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 ZHEEVD 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 241 of file zhegvd.f.

243*
244* -- LAPACK driver routine --
245* -- LAPACK is a software package provided by Univ. of Tennessee, --
246* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
247*
248* .. Scalar Arguments ..
249 CHARACTER JOBZ, UPLO
250 INTEGER INFO, ITYPE, LDA, LDB, LIWORK, LRWORK, LWORK, N
251* ..
252* .. Array Arguments ..
253 INTEGER IWORK( * )
254 DOUBLE PRECISION RWORK( * ), W( * )
255 COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( * )
256* ..
257*
258* =====================================================================
259*
260* .. Parameters ..
261 COMPLEX*16 CONE
262 parameter( cone = ( 1.0d+0, 0.0d+0 ) )
263* ..
264* .. Local Scalars ..
265 LOGICAL LQUERY, UPPER, WANTZ
266 CHARACTER TRANS
267 INTEGER LIOPT, LIWMIN, LOPT, LROPT, LRWMIN, LWMIN
268* ..
269* .. External Functions ..
270 LOGICAL LSAME
271 EXTERNAL lsame
272* ..
273* .. External Subroutines ..
274 EXTERNAL xerbla, zheevd, zhegst, zpotrf, ztrmm, ztrsm
275* ..
276* .. Intrinsic Functions ..
277 INTRINSIC dble, max
278* ..
279* .. Executable Statements ..
280*
281* Test the input parameters.
282*
283 wantz = lsame( jobz, 'V' )
284 upper = lsame( uplo, 'U' )
285 lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
286*
287 info = 0
288 IF( n.LE.1 ) THEN
289 lwmin = 1
290 lrwmin = 1
291 liwmin = 1
292 ELSE IF( wantz ) THEN
293 lwmin = 2*n + n*n
294 lrwmin = 1 + 5*n + 2*n*n
295 liwmin = 3 + 5*n
296 ELSE
297 lwmin = n + 1
298 lrwmin = n
299 liwmin = 1
300 END IF
301 lopt = lwmin
302 lropt = lrwmin
303 liopt = liwmin
304 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
305 info = -1
306 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
307 info = -2
308 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
309 info = -3
310 ELSE IF( n.LT.0 ) THEN
311 info = -4
312 ELSE IF( lda.LT.max( 1, n ) ) THEN
313 info = -6
314 ELSE IF( ldb.LT.max( 1, n ) ) THEN
315 info = -8
316 END IF
317*
318 IF( info.EQ.0 ) THEN
319 work( 1 ) = lopt
320 rwork( 1 ) = lropt
321 iwork( 1 ) = liopt
322*
323 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
324 info = -11
325 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
326 info = -13
327 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
328 info = -15
329 END IF
330 END IF
331*
332 IF( info.NE.0 ) THEN
333 CALL xerbla( 'ZHEGVD', -info )
334 RETURN
335 ELSE IF( lquery ) THEN
336 RETURN
337 END IF
338*
339* Quick return if possible
340*
341 IF( n.EQ.0 )
342 $ RETURN
343*
344* Form a Cholesky factorization of B.
345*
346 CALL zpotrf( uplo, n, b, ldb, info )
347 IF( info.NE.0 ) THEN
348 info = n + info
349 RETURN
350 END IF
351*
352* Transform problem to standard eigenvalue problem and solve.
353*
354 CALL zhegst( itype, uplo, n, a, lda, b, ldb, info )
355 CALL zheevd( jobz, uplo, n, a, lda, w, work, lwork, rwork, lrwork,
356 $ iwork, liwork, info )
357 lopt = int( max( dble( lopt ), dble( work( 1 ) ) ) )
358 lropt = int( max( dble( lropt ), dble( rwork( 1 ) ) ) )
359 liopt = int( max( dble( liopt ), dble( iwork( 1 ) ) ) )
360*
361 IF( wantz .AND. info.EQ.0 ) THEN
362*
363* Backtransform eigenvectors to the original problem.
364*
365 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
366*
367* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
368* backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
369*
370 IF( upper ) THEN
371 trans = 'N'
372 ELSE
373 trans = 'C'
374 END IF
375*
376 CALL ztrsm( 'Left', uplo, trans, 'Non-unit', n, n, cone,
377 $ b, ldb, a, lda )
378*
379 ELSE IF( itype.EQ.3 ) THEN
380*
381* For B*A*x=(lambda)*x;
382* backtransform eigenvectors: x = L*y or U**H *y
383*
384 IF( upper ) THEN
385 trans = 'C'
386 ELSE
387 trans = 'N'
388 END IF
389*
390 CALL ztrmm( 'Left', uplo, trans, 'Non-unit', n, n, cone,
391 $ b, ldb, a, lda )
392 END IF
393 END IF
394*
395 work( 1 ) = lopt
396 rwork( 1 ) = lropt
397 iwork( 1 ) = liopt
398*
399 RETURN
400*
401* End of ZHEGVD
402*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zheevd(jobz, uplo, n, a, lda, w, work, lwork, rwork, lrwork, iwork, liwork, info)
ZHEEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices
Definition zheevd.f:199
subroutine zhegst(itype, uplo, n, a, lda, b, ldb, info)
ZHEGST
Definition zhegst.f:128
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zpotrf(uplo, n, a, lda, info)
ZPOTRF
Definition zpotrf.f:107
subroutine ztrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRMM
Definition ztrmm.f:177
subroutine ztrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRSM
Definition ztrsm.f:180
Here is the call graph for this function:
Here is the caller graph for this function: