![]() ![]() | ![]() |
Mise à jour de lapack vers la version 3.3.0.
1: SUBROUTINE ZHEEVX( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, 2: $ ABSTOL, M, W, Z, LDZ, WORK, LWORK, RWORK, 3: $ IWORK, IFAIL, INFO ) 4: * 5: * -- LAPACK driver routine (version 3.2) -- 6: * -- LAPACK is a software package provided by Univ. of Tennessee, -- 7: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 8: * November 2006 9: * 10: * .. Scalar Arguments .. 11: CHARACTER JOBZ, RANGE, UPLO 12: INTEGER IL, INFO, IU, LDA, LDZ, LWORK, M, N 13: DOUBLE PRECISION ABSTOL, VL, VU 14: * .. 15: * .. Array Arguments .. 16: INTEGER IFAIL( * ), IWORK( * ) 17: DOUBLE PRECISION RWORK( * ), W( * ) 18: COMPLEX*16 A( LDA, * ), WORK( * ), Z( LDZ, * ) 19: * .. 20: * 21: * Purpose 22: * ======= 23: * 24: * ZHEEVX computes selected eigenvalues and, optionally, eigenvectors 25: * of a complex Hermitian matrix A. Eigenvalues and eigenvectors can 26: * be selected by specifying either a range of values or a range of 27: * indices for the desired eigenvalues. 28: * 29: * Arguments 30: * ========= 31: * 32: * JOBZ (input) CHARACTER*1 33: * = 'N': Compute eigenvalues only; 34: * = 'V': Compute eigenvalues and eigenvectors. 35: * 36: * RANGE (input) CHARACTER*1 37: * = 'A': all eigenvalues will be found. 38: * = 'V': all eigenvalues in the half-open interval (VL,VU] 39: * will be found. 40: * = 'I': the IL-th through IU-th eigenvalues will be found. 41: * 42: * UPLO (input) CHARACTER*1 43: * = 'U': Upper triangle of A is stored; 44: * = 'L': Lower triangle of A is stored. 45: * 46: * N (input) INTEGER 47: * The order of the matrix A. N >= 0. 48: * 49: * A (input/output) COMPLEX*16 array, dimension (LDA, N) 50: * On entry, the Hermitian matrix A. If UPLO = 'U', the 51: * leading N-by-N upper triangular part of A contains the 52: * upper triangular part of the matrix A. If UPLO = 'L', 53: * the leading N-by-N lower triangular part of A contains 54: * the lower triangular part of the matrix A. 55: * On exit, the lower triangle (if UPLO='L') or the upper 56: * triangle (if UPLO='U') of A, including the diagonal, is 57: * destroyed. 58: * 59: * LDA (input) INTEGER 60: * The leading dimension of the array A. LDA >= max(1,N). 61: * 62: * VL (input) DOUBLE PRECISION 63: * VU (input) DOUBLE PRECISION 64: * If RANGE='V', the lower and upper bounds of the interval to 65: * be searched for eigenvalues. VL < VU. 66: * Not referenced if RANGE = 'A' or 'I'. 67: * 68: * IL (input) INTEGER 69: * IU (input) INTEGER 70: * If RANGE='I', the indices (in ascending order) of the 71: * smallest and largest eigenvalues to be returned. 72: * 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. 73: * Not referenced if RANGE = 'A' or 'V'. 74: * 75: * ABSTOL (input) DOUBLE PRECISION 76: * The absolute error tolerance for the eigenvalues. 77: * An approximate eigenvalue is accepted as converged 78: * when it is determined to lie in an interval [a,b] 79: * of width less than or equal to 80: * 81: * ABSTOL + EPS * max( |a|,|b| ) , 82: * 83: * where EPS is the machine precision. If ABSTOL is less than 84: * or equal to zero, then EPS*|T| will be used in its place, 85: * where |T| is the 1-norm of the tridiagonal matrix obtained 86: * by reducing A to tridiagonal form. 87: * 88: * Eigenvalues will be computed most accurately when ABSTOL is 89: * set to twice the underflow threshold 2*DLAMCH('S'), not zero. 90: * If this routine returns with INFO>0, indicating that some 91: * eigenvectors did not converge, try setting ABSTOL to 92: * 2*DLAMCH('S'). 93: * 94: * See "Computing Small Singular Values of Bidiagonal Matrices 95: * with Guaranteed High Relative Accuracy," by Demmel and 96: * Kahan, LAPACK Working Note #3. 97: * 98: * M (output) INTEGER 99: * The total number of eigenvalues found. 0 <= M <= N. 100: * If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1. 101: * 102: * W (output) DOUBLE PRECISION array, dimension (N) 103: * On normal exit, the first M elements contain the selected 104: * eigenvalues in ascending order. 105: * 106: * Z (output) COMPLEX*16 array, dimension (LDZ, max(1,M)) 107: * If JOBZ = 'V', then if INFO = 0, the first M columns of Z 108: * contain the orthonormal eigenvectors of the matrix A 109: * corresponding to the selected eigenvalues, with the i-th 110: * column of Z holding the eigenvector associated with W(i). 111: * If an eigenvector fails to converge, then that column of Z 112: * contains the latest approximation to the eigenvector, and the 113: * index of the eigenvector is returned in IFAIL. 114: * If JOBZ = 'N', then Z is not referenced. 115: * Note: the user must ensure that at least max(1,M) columns are 116: * supplied in the array Z; if RANGE = 'V', the exact value of M 117: * is not known in advance and an upper bound must be used. 118: * 119: * LDZ (input) INTEGER 120: * The leading dimension of the array Z. LDZ >= 1, and if 121: * JOBZ = 'V', LDZ >= max(1,N). 122: * 123: * WORK (workspace/output) COMPLEX*16 array, dimension (MAX(1,LWORK)) 124: * On exit, if INFO = 0, WORK(1) returns the optimal LWORK. 125: * 126: * LWORK (input) INTEGER 127: * The length of the array WORK. LWORK >= 1, when N <= 1; 128: * otherwise 2*N. 129: * For optimal efficiency, LWORK >= (NB+1)*N, 130: * where NB is the max of the blocksize for ZHETRD and for 131: * ZUNMTR as returned by ILAENV. 132: * 133: * If LWORK = -1, then a workspace query is assumed; the routine 134: * only calculates the optimal size of the WORK array, returns 135: * this value as the first entry of the WORK array, and no error 136: * message related to LWORK is issued by XERBLA. 137: * 138: * RWORK (workspace) DOUBLE PRECISION array, dimension (7*N) 139: * 140: * IWORK (workspace) INTEGER array, dimension (5*N) 141: * 142: * IFAIL (output) INTEGER array, dimension (N) 143: * If JOBZ = 'V', then if INFO = 0, the first M elements of 144: * IFAIL are zero. If INFO > 0, then IFAIL contains the 145: * indices of the eigenvectors that failed to converge. 146: * If JOBZ = 'N', then IFAIL is not referenced. 147: * 148: * INFO (output) INTEGER 149: * = 0: successful exit 150: * < 0: if INFO = -i, the i-th argument had an illegal value 151: * > 0: if INFO = i, then i eigenvectors failed to converge. 152: * Their indices are stored in array IFAIL. 153: * 154: * ===================================================================== 155: * 156: * .. Parameters .. 157: DOUBLE PRECISION ZERO, ONE 158: PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 ) 159: COMPLEX*16 CONE 160: PARAMETER ( CONE = ( 1.0D+0, 0.0D+0 ) ) 161: * .. 162: * .. Local Scalars .. 163: LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, TEST, VALEIG, 164: $ WANTZ 165: CHARACTER ORDER 166: INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL, 167: $ INDISP, INDIWK, INDRWK, INDTAU, INDWRK, ISCALE, 168: $ ITMP1, J, JJ, LLWORK, LWKMIN, LWKOPT, NB, 169: $ NSPLIT 170: DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, 171: $ SIGMA, SMLNUM, TMP1, VLL, VUU 172: * .. 173: * .. External Functions .. 174: LOGICAL LSAME 175: INTEGER ILAENV 176: DOUBLE PRECISION DLAMCH, ZLANHE 177: EXTERNAL LSAME, ILAENV, DLAMCH, ZLANHE 178: * .. 179: * .. External Subroutines .. 180: EXTERNAL DCOPY, DSCAL, DSTEBZ, DSTERF, XERBLA, ZDSCAL, 181: $ ZHETRD, ZLACPY, ZSTEIN, ZSTEQR, ZSWAP, ZUNGTR, 182: $ ZUNMTR 183: * .. 184: * .. Intrinsic Functions .. 185: INTRINSIC DBLE, MAX, MIN, SQRT 186: * .. 187: * .. Executable Statements .. 188: * 189: * Test the input parameters. 190: * 191: LOWER = LSAME( UPLO, 'L' ) 192: WANTZ = LSAME( JOBZ, 'V' ) 193: ALLEIG = LSAME( RANGE, 'A' ) 194: VALEIG = LSAME( RANGE, 'V' ) 195: INDEIG = LSAME( RANGE, 'I' ) 196: LQUERY = ( LWORK.EQ.-1 ) 197: * 198: INFO = 0 199: IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN 200: INFO = -1 201: ELSE IF( .NOT.( ALLEIG .OR. VALEIG .OR. INDEIG ) ) THEN 202: INFO = -2 203: ELSE IF( .NOT.( LOWER .OR. LSAME( UPLO, 'U' ) ) ) THEN 204: INFO = -3 205: ELSE IF( N.LT.0 ) THEN 206: INFO = -4 207: ELSE IF( LDA.LT.MAX( 1, N ) ) THEN 208: INFO = -6 209: ELSE 210: IF( VALEIG ) THEN 211: IF( N.GT.0 .AND. VU.LE.VL ) 212: $ INFO = -8 213: ELSE IF( INDEIG ) THEN 214: IF( IL.LT.1 .OR. IL.GT.MAX( 1, N ) ) THEN 215: INFO = -9 216: ELSE IF( IU.LT.MIN( N, IL ) .OR. IU.GT.N ) THEN 217: INFO = -10 218: END IF 219: END IF 220: END IF 221: IF( INFO.EQ.0 ) THEN 222: IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN 223: INFO = -15 224: END IF 225: END IF 226: * 227: IF( INFO.EQ.0 ) THEN 228: IF( N.LE.1 ) THEN 229: LWKMIN = 1 230: WORK( 1 ) = LWKMIN 231: ELSE 232: LWKMIN = 2*N 233: NB = ILAENV( 1, 'ZHETRD', UPLO, N, -1, -1, -1 ) 234: NB = MAX( NB, ILAENV( 1, 'ZUNMTR', UPLO, N, -1, -1, -1 ) ) 235: LWKOPT = MAX( 1, ( NB + 1 )*N ) 236: WORK( 1 ) = LWKOPT 237: END IF 238: * 239: IF( LWORK.LT.MAX( 1, 2*N ) .AND. .NOT.LQUERY ) 240: $ INFO = -17 241: END IF 242: * 243: IF( INFO.NE.0 ) THEN 244: CALL XERBLA( 'ZHEEVX', -INFO ) 245: RETURN 246: ELSE IF( LQUERY ) THEN 247: RETURN 248: END IF 249: * 250: * Quick return if possible 251: * 252: M = 0 253: IF( N.EQ.0 ) THEN 254: RETURN 255: END IF 256: * 257: IF( N.EQ.1 ) THEN 258: IF( ALLEIG .OR. INDEIG ) THEN 259: M = 1 260: W( 1 ) = A( 1, 1 ) 261: ELSE IF( VALEIG ) THEN 262: IF( VL.LT.DBLE( A( 1, 1 ) ) .AND. VU.GE.DBLE( A( 1, 1 ) ) ) 263: $ THEN 264: M = 1 265: W( 1 ) = A( 1, 1 ) 266: END IF 267: END IF 268: IF( WANTZ ) 269: $ Z( 1, 1 ) = CONE 270: RETURN 271: END IF 272: * 273: * Get machine constants. 274: * 275: SAFMIN = DLAMCH( 'Safe minimum' ) 276: EPS = DLAMCH( 'Precision' ) 277: SMLNUM = SAFMIN / EPS 278: BIGNUM = ONE / SMLNUM 279: RMIN = SQRT( SMLNUM ) 280: RMAX = MIN( SQRT( BIGNUM ), ONE / SQRT( SQRT( SAFMIN ) ) ) 281: * 282: * Scale matrix to allowable range, if necessary. 283: * 284: ISCALE = 0 285: ABSTLL = ABSTOL 286: IF( VALEIG ) THEN 287: VLL = VL 288: VUU = VU 289: END IF 290: ANRM = ZLANHE( 'M', UPLO, N, A, LDA, RWORK ) 291: IF( ANRM.GT.ZERO .AND. ANRM.LT.RMIN ) THEN 292: ISCALE = 1 293: SIGMA = RMIN / ANRM 294: ELSE IF( ANRM.GT.RMAX ) THEN 295: ISCALE = 1 296: SIGMA = RMAX / ANRM 297: END IF 298: IF( ISCALE.EQ.1 ) THEN 299: IF( LOWER ) THEN 300: DO 10 J = 1, N 301: CALL ZDSCAL( N-J+1, SIGMA, A( J, J ), 1 ) 302: 10 CONTINUE 303: ELSE 304: DO 20 J = 1, N 305: CALL ZDSCAL( J, SIGMA, A( 1, J ), 1 ) 306: 20 CONTINUE 307: END IF 308: IF( ABSTOL.GT.0 ) 309: $ ABSTLL = ABSTOL*SIGMA 310: IF( VALEIG ) THEN 311: VLL = VL*SIGMA 312: VUU = VU*SIGMA 313: END IF 314: END IF 315: * 316: * Call ZHETRD to reduce Hermitian matrix to tridiagonal form. 317: * 318: INDD = 1 319: INDE = INDD + N 320: INDRWK = INDE + N 321: INDTAU = 1 322: INDWRK = INDTAU + N 323: LLWORK = LWORK - INDWRK + 1 324: CALL ZHETRD( UPLO, N, A, LDA, RWORK( INDD ), RWORK( INDE ), 325: $ WORK( INDTAU ), WORK( INDWRK ), LLWORK, IINFO ) 326: * 327: * If all eigenvalues are desired and ABSTOL is less than or equal to 328: * zero, then call DSTERF or ZUNGTR and ZSTEQR. If this fails for 329: * some eigenvalue, then try DSTEBZ. 330: * 331: TEST = .FALSE. 332: IF( INDEIG ) THEN 333: IF( IL.EQ.1 .AND. IU.EQ.N ) THEN 334: TEST = .TRUE. 335: END IF 336: END IF 337: IF( ( ALLEIG .OR. TEST ) .AND. ( ABSTOL.LE.ZERO ) ) THEN 338: CALL DCOPY( N, RWORK( INDD ), 1, W, 1 ) 339: INDEE = INDRWK + 2*N 340: IF( .NOT.WANTZ ) THEN 341: CALL DCOPY( N-1, RWORK( INDE ), 1, RWORK( INDEE ), 1 ) 342: CALL DSTERF( N, W, RWORK( INDEE ), INFO ) 343: ELSE 344: CALL ZLACPY( 'A', N, N, A, LDA, Z, LDZ ) 345: CALL ZUNGTR( UPLO, N, Z, LDZ, WORK( INDTAU ), 346: $ WORK( INDWRK ), LLWORK, IINFO ) 347: CALL DCOPY( N-1, RWORK( INDE ), 1, RWORK( INDEE ), 1 ) 348: CALL ZSTEQR( JOBZ, N, W, RWORK( INDEE ), Z, LDZ, 349: $ RWORK( INDRWK ), INFO ) 350: IF( INFO.EQ.0 ) THEN 351: DO 30 I = 1, N 352: IFAIL( I ) = 0 353: 30 CONTINUE 354: END IF 355: END IF 356: IF( INFO.EQ.0 ) THEN 357: M = N 358: GO TO 40 359: END IF 360: INFO = 0 361: END IF 362: * 363: * Otherwise, call DSTEBZ and, if eigenvectors are desired, ZSTEIN. 364: * 365: IF( WANTZ ) THEN 366: ORDER = 'B' 367: ELSE 368: ORDER = 'E' 369: END IF 370: INDIBL = 1 371: INDISP = INDIBL + N 372: INDIWK = INDISP + N 373: CALL DSTEBZ( RANGE, ORDER, N, VLL, VUU, IL, IU, ABSTLL, 374: $ RWORK( INDD ), RWORK( INDE ), M, NSPLIT, W, 375: $ IWORK( INDIBL ), IWORK( INDISP ), RWORK( INDRWK ), 376: $ IWORK( INDIWK ), INFO ) 377: * 378: IF( WANTZ ) THEN 379: CALL ZSTEIN( N, RWORK( INDD ), RWORK( INDE ), M, W, 380: $ IWORK( INDIBL ), IWORK( INDISP ), Z, LDZ, 381: $ RWORK( INDRWK ), IWORK( INDIWK ), IFAIL, INFO ) 382: * 383: * Apply unitary matrix used in reduction to tridiagonal 384: * form to eigenvectors returned by ZSTEIN. 385: * 386: CALL ZUNMTR( 'L', UPLO, 'N', N, M, A, LDA, WORK( INDTAU ), Z, 387: $ LDZ, WORK( INDWRK ), LLWORK, IINFO ) 388: END IF 389: * 390: * If matrix was scaled, then rescale eigenvalues appropriately. 391: * 392: 40 CONTINUE 393: IF( ISCALE.EQ.1 ) THEN 394: IF( INFO.EQ.0 ) THEN 395: IMAX = M 396: ELSE 397: IMAX = INFO - 1 398: END IF 399: CALL DSCAL( IMAX, ONE / SIGMA, W, 1 ) 400: END IF 401: * 402: * If eigenvalues are not in order, then sort them, along with 403: * eigenvectors. 404: * 405: IF( WANTZ ) THEN 406: DO 60 J = 1, M - 1 407: I = 0 408: TMP1 = W( J ) 409: DO 50 JJ = J + 1, M 410: IF( W( JJ ).LT.TMP1 ) THEN 411: I = JJ 412: TMP1 = W( JJ ) 413: END IF 414: 50 CONTINUE 415: * 416: IF( I.NE.0 ) THEN 417: ITMP1 = IWORK( INDIBL+I-1 ) 418: W( I ) = W( J ) 419: IWORK( INDIBL+I-1 ) = IWORK( INDIBL+J-1 ) 420: W( J ) = TMP1 421: IWORK( INDIBL+J-1 ) = ITMP1 422: CALL ZSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 ) 423: IF( INFO.NE.0 ) THEN 424: ITMP1 = IFAIL( I ) 425: IFAIL( I ) = IFAIL( J ) 426: IFAIL( J ) = ITMP1 427: END IF 428: END IF 429: 60 CONTINUE 430: END IF 431: * 432: * Set WORK(1) to optimal complex workspace size. 433: * 434: WORK( 1 ) = LWKOPT 435: * 436: RETURN 437: * 438: * End of ZHEEVX 439: * 440: END