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

◆ zhpgvd()

subroutine zhpgvd ( integer  ITYPE,
character  JOBZ,
character  UPLO,
integer  N,
complex*16, dimension( * )  AP,
complex*16, dimension( * )  BP,
double precision, dimension( * )  W,
complex*16, dimension( ldz, * )  Z,
integer  LDZ,
complex*16, dimension( * )  WORK,
integer  LWORK,
double precision, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

ZHPGVD

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

Purpose:
 ZHPGVD 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, stored in packed format, 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]AP
          AP is COMPLEX*16 array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the Hermitian matrix
          A, packed columnwise in a linear array.  The j-th column of A
          is stored in the array AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.

          On exit, the contents of AP are destroyed.
[in,out]BP
          BP is COMPLEX*16 array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the Hermitian matrix
          B, packed columnwise in a linear array.  The j-th column of B
          is stored in the array BP as follows:
          if UPLO = 'U', BP(i + (j-1)*j/2) = B(i,j) for 1<=i<=j;
          if UPLO = 'L', BP(i + (j-1)*(2*n-j)/2) = B(i,j) for j<=i<=n.

          On exit, the triangular factor U or L from the Cholesky
          factorization B = U**H*U or B = L*L**H, in the same storage
          format as B.
[out]W
          W is DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is COMPLEX*16 array, dimension (LDZ, N)
          If JOBZ = 'V', then if INFO = 0, Z 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 Z is not referenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= 1, and if
          JOBZ = 'V', LDZ >= max(1,N).
[out]WORK
          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the required LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.
          If N <= 1,               LWORK >= 1.
          If JOBZ = 'N' and N > 1, LWORK >= N.
          If JOBZ = 'V' and N > 1, LWORK >= 2*N.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the required 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 required LRWORK.
[in]LRWORK
          LRWORK is INTEGER
          The dimension of 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 required 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 required LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of array IWORK.
          If JOBZ  = 'N' or 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 required 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:  ZPPTRF or ZHPEVD returned an error code:
             <= N:  if INFO = i, ZHPEVD failed to converge;
                    i off-diagonal elements of an intermediate
                    tridiagonal form did not convergeto 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.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 229 of file zhpgvd.f.

231*
232* -- LAPACK driver routine --
233* -- LAPACK is a software package provided by Univ. of Tennessee, --
234* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
235*
236* .. Scalar Arguments ..
237 CHARACTER JOBZ, UPLO
238 INTEGER INFO, ITYPE, LDZ, LIWORK, LRWORK, LWORK, N
239* ..
240* .. Array Arguments ..
241 INTEGER IWORK( * )
242 DOUBLE PRECISION RWORK( * ), W( * )
243 COMPLEX*16 AP( * ), BP( * ), WORK( * ), Z( LDZ, * )
244* ..
245*
246* =====================================================================
247*
248* .. Local Scalars ..
249 LOGICAL LQUERY, UPPER, WANTZ
250 CHARACTER TRANS
251 INTEGER J, LIWMIN, LRWMIN, LWMIN, NEIG
252* ..
253* .. External Functions ..
254 LOGICAL LSAME
255 EXTERNAL lsame
256* ..
257* .. External Subroutines ..
258 EXTERNAL xerbla, zhpevd, zhpgst, zpptrf, ztpmv, ztpsv
259* ..
260* .. Intrinsic Functions ..
261 INTRINSIC dble, max
262* ..
263* .. Executable Statements ..
264*
265* Test the input parameters.
266*
267 wantz = lsame( jobz, 'V' )
268 upper = lsame( uplo, 'U' )
269 lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
270*
271 info = 0
272 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
273 info = -1
274 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
275 info = -2
276 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
277 info = -3
278 ELSE IF( n.LT.0 ) THEN
279 info = -4
280 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
281 info = -9
282 END IF
283*
284 IF( info.EQ.0 ) THEN
285 IF( n.LE.1 ) THEN
286 lwmin = 1
287 liwmin = 1
288 lrwmin = 1
289 ELSE
290 IF( wantz ) THEN
291 lwmin = 2*n
292 lrwmin = 1 + 5*n + 2*n**2
293 liwmin = 3 + 5*n
294 ELSE
295 lwmin = n
296 lrwmin = n
297 liwmin = 1
298 END IF
299 END IF
300*
301 work( 1 ) = lwmin
302 rwork( 1 ) = lrwmin
303 iwork( 1 ) = liwmin
304 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
305 info = -11
306 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
307 info = -13
308 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
309 info = -15
310 END IF
311 END IF
312*
313 IF( info.NE.0 ) THEN
314 CALL xerbla( 'ZHPGVD', -info )
315 RETURN
316 ELSE IF( lquery ) THEN
317 RETURN
318 END IF
319*
320* Quick return if possible
321*
322 IF( n.EQ.0 )
323 $ RETURN
324*
325* Form a Cholesky factorization of B.
326*
327 CALL zpptrf( uplo, n, bp, info )
328 IF( info.NE.0 ) THEN
329 info = n + info
330 RETURN
331 END IF
332*
333* Transform problem to standard eigenvalue problem and solve.
334*
335 CALL zhpgst( itype, uplo, n, ap, bp, info )
336 CALL zhpevd( jobz, uplo, n, ap, w, z, ldz, work, lwork, rwork,
337 $ lrwork, iwork, liwork, info )
338 lwmin = int( max( dble( lwmin ), dble( work( 1 ) ) ) )
339 lrwmin = int( max( dble( lrwmin ), dble( rwork( 1 ) ) ) )
340 liwmin = int( max( dble( liwmin ), dble( iwork( 1 ) ) ) )
341*
342 IF( wantz ) THEN
343*
344* Backtransform eigenvectors to the original problem.
345*
346 neig = n
347 IF( info.GT.0 )
348 $ neig = info - 1
349 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
350*
351* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
352* backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
353*
354 IF( upper ) THEN
355 trans = 'N'
356 ELSE
357 trans = 'C'
358 END IF
359*
360 DO 10 j = 1, neig
361 CALL ztpsv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
362 $ 1 )
363 10 CONTINUE
364*
365 ELSE IF( itype.EQ.3 ) THEN
366*
367* For B*A*x=(lambda)*x;
368* backtransform eigenvectors: x = L*y or U**H *y
369*
370 IF( upper ) THEN
371 trans = 'C'
372 ELSE
373 trans = 'N'
374 END IF
375*
376 DO 20 j = 1, neig
377 CALL ztpmv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
378 $ 1 )
379 20 CONTINUE
380 END IF
381 END IF
382*
383 work( 1 ) = lwmin
384 rwork( 1 ) = lrwmin
385 iwork( 1 ) = liwmin
386 RETURN
387*
388* End of ZHPGVD
389*
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine ztpsv(UPLO, TRANS, DIAG, N, AP, X, INCX)
ZTPSV
Definition: ztpsv.f:144
subroutine ztpmv(UPLO, TRANS, DIAG, N, AP, X, INCX)
ZTPMV
Definition: ztpmv.f:142
subroutine zhpgst(ITYPE, UPLO, N, AP, BP, INFO)
ZHPGST
Definition: zhpgst.f:113
subroutine zpptrf(UPLO, N, AP, INFO)
ZPPTRF
Definition: zpptrf.f:119
subroutine zhpevd(JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
ZHPEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: zhpevd.f:200
Here is the call graph for this function:
Here is the caller graph for this function: