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

◆ chegvd()

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.
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
                    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 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 241 of file chegvd.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 REAL RWORK( * ), W( * )
255 COMPLEX A( LDA, * ), B( LDB, * ), WORK( * )
256* ..
257*
258* =====================================================================
259*
260* .. Parameters ..
261 COMPLEX CONE
262 parameter( cone = ( 1.0e+0, 0.0e+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 REAL SROUNDUP_LWORK
272 EXTERNAL lsame, sroundup_lwork
273* ..
274* .. External Subroutines ..
275 EXTERNAL cheevd, chegst, cpotrf, ctrmm, ctrsm, xerbla
276* ..
277* .. Intrinsic Functions ..
278 INTRINSIC max, real
279* ..
280* .. Executable Statements ..
281*
282* Test the input parameters.
283*
284 wantz = lsame( jobz, 'V' )
285 upper = lsame( uplo, 'U' )
286 lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
287*
288 info = 0
289 IF( n.LE.1 ) THEN
290 lwmin = 1
291 lrwmin = 1
292 liwmin = 1
293 ELSE IF( wantz ) THEN
294 lwmin = 2*n + n*n
295 lrwmin = 1 + 5*n + 2*n*n
296 liwmin = 3 + 5*n
297 ELSE
298 lwmin = n + 1
299 lrwmin = n
300 liwmin = 1
301 END IF
302 lopt = lwmin
303 lropt = lrwmin
304 liopt = liwmin
305 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
306 info = -1
307 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
308 info = -2
309 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
310 info = -3
311 ELSE IF( n.LT.0 ) THEN
312 info = -4
313 ELSE IF( lda.LT.max( 1, n ) ) THEN
314 info = -6
315 ELSE IF( ldb.LT.max( 1, n ) ) THEN
316 info = -8
317 END IF
318*
319 IF( info.EQ.0 ) THEN
320 work( 1 ) = sroundup_lwork(lopt)
321 rwork( 1 ) = lropt
322 iwork( 1 ) = liopt
323*
324 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
325 info = -11
326 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
327 info = -13
328 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
329 info = -15
330 END IF
331 END IF
332*
333 IF( info.NE.0 ) THEN
334 CALL xerbla( 'CHEGVD', -info )
335 RETURN
336 ELSE IF( lquery ) THEN
337 RETURN
338 END IF
339*
340* Quick return if possible
341*
342 IF( n.EQ.0 )
343 $ RETURN
344*
345* Form a Cholesky factorization of B.
346*
347 CALL cpotrf( uplo, n, b, ldb, info )
348 IF( info.NE.0 ) THEN
349 info = n + info
350 RETURN
351 END IF
352*
353* Transform problem to standard eigenvalue problem and solve.
354*
355 CALL chegst( itype, uplo, n, a, lda, b, ldb, info )
356 CALL cheevd( jobz, uplo, n, a, lda, w, work, lwork, rwork, lrwork,
357 $ iwork, liwork, info )
358 lopt = int( max( real( lopt ), real( work( 1 ) ) ) )
359 lropt = int( max( real( lropt ), real( rwork( 1 ) ) ) )
360 liopt = int( max( real( liopt ), real( iwork( 1 ) ) ) )
361*
362 IF( wantz .AND. info.EQ.0 ) THEN
363*
364* Backtransform eigenvectors to the original problem.
365*
366 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
367*
368* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
369* backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
370*
371 IF( upper ) THEN
372 trans = 'N'
373 ELSE
374 trans = 'C'
375 END IF
376*
377 CALL ctrsm( 'Left', uplo, trans, 'Non-unit', n, n, cone,
378 $ b, ldb, a, lda )
379*
380 ELSE IF( itype.EQ.3 ) THEN
381*
382* For B*A*x=(lambda)*x;
383* backtransform eigenvectors: x = L*y or U**H *y
384*
385 IF( upper ) THEN
386 trans = 'C'
387 ELSE
388 trans = 'N'
389 END IF
390*
391 CALL ctrmm( 'Left', uplo, trans, 'Non-unit', n, n, cone,
392 $ b, ldb, a, lda )
393 END IF
394 END IF
395*
396 work( 1 ) = sroundup_lwork(lopt)
397 rwork( 1 ) = lropt
398 iwork( 1 ) = liopt
399*
400 RETURN
401*
402* End of CHEGVD
403*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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:199
subroutine chegst(itype, uplo, n, a, lda, b, ldb, info)
CHEGST
Definition chegst.f:128
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine cpotrf(uplo, n, a, lda, info)
CPOTRF
Definition cpotrf.f:107
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine ctrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
CTRMM
Definition ctrmm.f:177
subroutine ctrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
CTRSM
Definition ctrsm.f:180
Here is the call graph for this function:
Here is the caller graph for this function: