LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zhegv ( 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  INFO 
)

ZHEGV

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

Purpose:
 ZHEGV 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.
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 positive definite 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.  LWORK >= max(1,2*N-1).
          For optimal efficiency, LWORK >= (NB+1)*N,
          where NB is the blocksize for ZHETRD returned by ILAENV.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK is issued by XERBLA.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (max(1, 3*N-2))
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  ZPOTRF or ZHEEV returned an error code:
             <= N:  if INFO = i, ZHEEV failed to converge;
                    i off-diagonal elements of an intermediate
                    tridiagonal form did not converge to zero;
             > 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

Definition at line 183 of file zhegv.f.

183 *
184 * -- LAPACK driver routine (version 3.6.0) --
185 * -- LAPACK is a software package provided by Univ. of Tennessee, --
186 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
187 * November 2015
188 *
189 * .. Scalar Arguments ..
190  CHARACTER jobz, uplo
191  INTEGER info, itype, lda, ldb, lwork, n
192 * ..
193 * .. Array Arguments ..
194  DOUBLE PRECISION rwork( * ), w( * )
195  COMPLEX*16 a( lda, * ), b( ldb, * ), work( * )
196 * ..
197 *
198 * =====================================================================
199 *
200 * .. Parameters ..
201  COMPLEX*16 one
202  parameter ( one = ( 1.0d+0, 0.0d+0 ) )
203 * ..
204 * .. Local Scalars ..
205  LOGICAL lquery, upper, wantz
206  CHARACTER trans
207  INTEGER lwkopt, nb, neig
208 * ..
209 * .. External Functions ..
210  LOGICAL lsame
211  INTEGER ilaenv
212  EXTERNAL lsame, ilaenv
213 * ..
214 * .. External Subroutines ..
215  EXTERNAL xerbla, zheev, zhegst, zpotrf, ztrmm, ztrsm
216 * ..
217 * .. Intrinsic Functions ..
218  INTRINSIC max
219 * ..
220 * .. Executable Statements ..
221 *
222 * Test the input parameters.
223 *
224  wantz = lsame( jobz, 'V' )
225  upper = lsame( uplo, 'U' )
226  lquery = ( lwork.EQ.-1 )
227 *
228  info = 0
229  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
230  info = -1
231  ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
232  info = -2
233  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
234  info = -3
235  ELSE IF( n.LT.0 ) THEN
236  info = -4
237  ELSE IF( lda.LT.max( 1, n ) ) THEN
238  info = -6
239  ELSE IF( ldb.LT.max( 1, n ) ) THEN
240  info = -8
241  END IF
242 *
243  IF( info.EQ.0 ) THEN
244  nb = ilaenv( 1, 'ZHETRD', uplo, n, -1, -1, -1 )
245  lwkopt = max( 1, ( nb + 1 )*n )
246  work( 1 ) = lwkopt
247 *
248  IF( lwork.LT.max( 1, 2*n - 1 ) .AND. .NOT.lquery ) THEN
249  info = -11
250  END IF
251  END IF
252 *
253  IF( info.NE.0 ) THEN
254  CALL xerbla( 'ZHEGV ', -info )
255  RETURN
256  ELSE IF( lquery ) THEN
257  RETURN
258  END IF
259 *
260 * Quick return if possible
261 *
262  IF( n.EQ.0 )
263  $ RETURN
264 *
265 * Form a Cholesky factorization of B.
266 *
267  CALL zpotrf( uplo, n, b, ldb, info )
268  IF( info.NE.0 ) THEN
269  info = n + info
270  RETURN
271  END IF
272 *
273 * Transform problem to standard eigenvalue problem and solve.
274 *
275  CALL zhegst( itype, uplo, n, a, lda, b, ldb, info )
276  CALL zheev( jobz, uplo, n, a, lda, w, work, lwork, rwork, info )
277 *
278  IF( wantz ) THEN
279 *
280 * Backtransform eigenvectors to the original problem.
281 *
282  neig = n
283  IF( info.GT.0 )
284  $ neig = info - 1
285  IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
286 *
287 * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
288 * backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
289 *
290  IF( upper ) THEN
291  trans = 'N'
292  ELSE
293  trans = 'C'
294  END IF
295 *
296  CALL ztrsm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
297  $ b, ldb, a, lda )
298 *
299  ELSE IF( itype.EQ.3 ) THEN
300 *
301 * For B*A*x=(lambda)*x;
302 * backtransform eigenvectors: x = L*y or U**H *y
303 *
304  IF( upper ) THEN
305  trans = 'C'
306  ELSE
307  trans = 'N'
308  END IF
309 *
310  CALL ztrmm( 'Left', uplo, trans, 'Non-unit', n, neig, one,
311  $ b, ldb, a, lda )
312  END IF
313  END IF
314 *
315  work( 1 ) = lwkopt
316 *
317  RETURN
318 *
319 * End of ZHEGV
320 *
subroutine zpotrf(UPLO, N, A, LDA, INFO)
ZPOTRF VARIANT: right looking block version of the algorithm, calling Level 3 BLAS.
Definition: zpotrf.f:102
subroutine zhegst(ITYPE, UPLO, N, A, LDA, B, LDB, INFO)
ZHEGST
Definition: zhegst.f:129
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ztrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
ZTRMM
Definition: ztrmm.f:179
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
subroutine zheev(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, RWORK, INFO)
ZHEEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices ...
Definition: zheev.f:142
subroutine ztrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
ZTRSM
Definition: ztrsm.f:182
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: