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

◆ zhpgvx()

subroutine zhpgvx ( integer  itype,
character  jobz,
character  range,
character  uplo,
integer  n,
complex*16, dimension( * )  ap,
complex*16, dimension( * )  bp,
double precision  vl,
double precision  vu,
integer  il,
integer  iu,
double precision  abstol,
integer  m,
double precision, dimension( * )  w,
complex*16, dimension( ldz, * )  z,
integer  ldz,
complex*16, dimension( * )  work,
double precision, dimension( * )  rwork,
integer, dimension( * )  iwork,
integer, dimension( * )  ifail,
integer  info 
)

ZHPGVX

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

Purpose:
 ZHPGVX computes selected eigenvalues and, optionally, 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.  Eigenvalues and eigenvectors can be selected by
 specifying either a range of values or a range of indices for the
 desired eigenvalues.
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]RANGE
          RANGE is CHARACTER*1
          = 'A': all eigenvalues will be found;
          = 'V': all eigenvalues in the half-open interval (VL,VU]
                 will be found;
          = 'I': the IL-th through IU-th eigenvalues will be found.
[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.
[in]VL
          VL is DOUBLE PRECISION

          If RANGE='V', the lower bound of the interval to
          be searched for eigenvalues. VL < VU.
          Not referenced if RANGE = 'A' or 'I'.
[in]VU
          VU is DOUBLE PRECISION

          If RANGE='V', the upper bound of the interval to
          be searched for eigenvalues. VL < VU.
          Not referenced if RANGE = 'A' or 'I'.
[in]IL
          IL is INTEGER

          If RANGE='I', the index of the
          smallest eigenvalue to be returned.
          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
          Not referenced if RANGE = 'A' or 'V'.
[in]IU
          IU is INTEGER

          If RANGE='I', the index of the
          largest eigenvalue to be returned.
          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
          Not referenced if RANGE = 'A' or 'V'.
[in]ABSTOL
          ABSTOL is DOUBLE PRECISION
          The absolute error tolerance for the eigenvalues.
          An approximate eigenvalue is accepted as converged
          when it is determined to lie in an interval [a,b]
          of width less than or equal to

                  ABSTOL + EPS *   max( |a|,|b| ) ,

          where EPS is the machine precision.  If ABSTOL is less than
          or equal to zero, then  EPS*|T|  will be used in its place,
          where |T| is the 1-norm of the tridiagonal matrix obtained
          by reducing AP to tridiagonal form.

          Eigenvalues will be computed most accurately when ABSTOL is
          set to twice the underflow threshold 2*DLAMCH('S'), not zero.
          If this routine returns with INFO>0, indicating that some
          eigenvectors did not converge, try setting ABSTOL to
          2*DLAMCH('S').
[out]M
          M is INTEGER
          The total number of eigenvalues found.  0 <= M <= N.
          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
[out]W
          W is DOUBLE PRECISION array, dimension (N)
          On normal exit, the first M elements contain the selected
          eigenvalues in ascending order.
[out]Z
          Z is COMPLEX*16 array, dimension (LDZ, N)
          If JOBZ = 'N', then Z is not referenced.
          If JOBZ = 'V', then if INFO = 0, the first M columns of Z
          contain the orthonormal eigenvectors of the matrix A
          corresponding to the selected eigenvalues, with the i-th
          column of Z holding the eigenvector associated with W(i).
          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 an eigenvector fails to converge, then that column of Z
          contains the latest approximation to the eigenvector, and the
          index of the eigenvector is returned in IFAIL.
          Note: the user must ensure that at least max(1,M) columns are
          supplied in the array Z; if RANGE = 'V', the exact value of M
          is not known in advance and an upper bound must be used.
[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 (2*N)
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (7*N)
[out]IWORK
          IWORK is INTEGER array, dimension (5*N)
[out]IFAIL
          IFAIL is INTEGER array, dimension (N)
          If JOBZ = 'V', then if INFO = 0, the first M elements of
          IFAIL are zero.  If INFO > 0, then IFAIL contains the
          indices of the eigenvectors that failed to converge.
          If JOBZ = 'N', then IFAIL is not referenced.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  ZPPTRF or ZHPEVX returned an error code:
             <= N:  if INFO = i, ZHPEVX failed to converge;
                    i eigenvectors failed to converge.  Their indices
                    are stored in array IFAIL.
             > 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.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 274 of file zhpgvx.f.

277*
278* -- LAPACK driver routine --
279* -- LAPACK is a software package provided by Univ. of Tennessee, --
280* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
281*
282* .. Scalar Arguments ..
283 CHARACTER JOBZ, RANGE, UPLO
284 INTEGER IL, INFO, ITYPE, IU, LDZ, M, N
285 DOUBLE PRECISION ABSTOL, VL, VU
286* ..
287* .. Array Arguments ..
288 INTEGER IFAIL( * ), IWORK( * )
289 DOUBLE PRECISION RWORK( * ), W( * )
290 COMPLEX*16 AP( * ), BP( * ), WORK( * ), Z( LDZ, * )
291* ..
292*
293* =====================================================================
294*
295* .. Local Scalars ..
296 LOGICAL ALLEIG, INDEIG, UPPER, VALEIG, WANTZ
297 CHARACTER TRANS
298 INTEGER J
299* ..
300* .. External Functions ..
301 LOGICAL LSAME
302 EXTERNAL lsame
303* ..
304* .. External Subroutines ..
305 EXTERNAL xerbla, zhpevx, zhpgst, zpptrf, ztpmv, ztpsv
306* ..
307* .. Intrinsic Functions ..
308 INTRINSIC min
309* ..
310* .. Executable Statements ..
311*
312* Test the input parameters.
313*
314 wantz = lsame( jobz, 'V' )
315 upper = lsame( uplo, 'U' )
316 alleig = lsame( range, 'A' )
317 valeig = lsame( range, 'V' )
318 indeig = lsame( range, 'I' )
319*
320 info = 0
321 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
322 info = -1
323 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
324 info = -2
325 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
326 info = -3
327 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
328 info = -4
329 ELSE IF( n.LT.0 ) THEN
330 info = -5
331 ELSE
332 IF( valeig ) THEN
333 IF( n.GT.0 .AND. vu.LE.vl ) THEN
334 info = -9
335 END IF
336 ELSE IF( indeig ) THEN
337 IF( il.LT.1 ) THEN
338 info = -10
339 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
340 info = -11
341 END IF
342 END IF
343 END IF
344 IF( info.EQ.0 ) THEN
345 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
346 info = -16
347 END IF
348 END IF
349*
350 IF( info.NE.0 ) THEN
351 CALL xerbla( 'ZHPGVX', -info )
352 RETURN
353 END IF
354*
355* Quick return if possible
356*
357 IF( n.EQ.0 )
358 $ RETURN
359*
360* Form a Cholesky factorization of B.
361*
362 CALL zpptrf( uplo, n, bp, info )
363 IF( info.NE.0 ) THEN
364 info = n + info
365 RETURN
366 END IF
367*
368* Transform problem to standard eigenvalue problem and solve.
369*
370 CALL zhpgst( itype, uplo, n, ap, bp, info )
371 CALL zhpevx( jobz, range, uplo, n, ap, vl, vu, il, iu, abstol, m,
372 $ w, z, ldz, work, rwork, iwork, ifail, info )
373*
374 IF( wantz ) THEN
375*
376* Backtransform eigenvectors to the original problem.
377*
378 IF( info.GT.0 )
379 $ m = info - 1
380 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
381*
382* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
383* backtransform eigenvectors: x = inv(L)**H *y or inv(U)*y
384*
385 IF( upper ) THEN
386 trans = 'N'
387 ELSE
388 trans = 'C'
389 END IF
390*
391 DO 10 j = 1, m
392 CALL ztpsv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
393 $ 1 )
394 10 CONTINUE
395*
396 ELSE IF( itype.EQ.3 ) THEN
397*
398* For B*A*x=(lambda)*x;
399* backtransform eigenvectors: x = L*y or U**H *y
400*
401 IF( upper ) THEN
402 trans = 'C'
403 ELSE
404 trans = 'N'
405 END IF
406*
407 DO 20 j = 1, m
408 CALL ztpmv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
409 $ 1 )
410 20 CONTINUE
411 END IF
412 END IF
413*
414 RETURN
415*
416* End of ZHPGVX
417*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zhpevx(jobz, range, uplo, n, ap, vl, vu, il, iu, abstol, m, w, z, ldz, work, rwork, iwork, ifail, info)
ZHPEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition zhpevx.f:240
subroutine zhpgst(itype, uplo, n, ap, bp, info)
ZHPGST
Definition zhpgst.f:113
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zpptrf(uplo, n, ap, info)
ZPPTRF
Definition zpptrf.f:119
subroutine ztpmv(uplo, trans, diag, n, ap, x, incx)
ZTPMV
Definition ztpmv.f:142
subroutine ztpsv(uplo, trans, diag, n, ap, x, incx)
ZTPSV
Definition ztpsv.f:144
Here is the call graph for this function:
Here is the caller graph for this function: