1*> \brief \b ZSTEMR 2* 3* =========== DOCUMENTATION =========== 4* 5* Online html documentation available at 6* http://www.netlib.org/lapack/explore-html/ 7* 8*> \htmlonly 9*> Download ZSTEMR + dependencies 10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zstemr.f"> 11*> [TGZ]</a> 12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zstemr.f"> 13*> [ZIP]</a> 14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zstemr.f"> 15*> [TXT]</a> 16*> \endhtmlonly 17* 18* Definition: 19* =========== 20* 21* SUBROUTINE ZSTEMR( JOBZ, RANGE, N, D, E, VL, VU, IL, IU, 22* M, W, Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, 23* IWORK, LIWORK, INFO ) 24* 25* .. Scalar Arguments .. 26* CHARACTER JOBZ, RANGE 27* LOGICAL TRYRAC 28* INTEGER IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N 29* DOUBLE PRECISION VL, VU 30* .. 31* .. Array Arguments .. 32* INTEGER ISUPPZ( * ), IWORK( * ) 33* DOUBLE PRECISION D( * ), E( * ), W( * ), WORK( * ) 34* COMPLEX*16 Z( LDZ, * ) 35* .. 36* 37* 38*> \par Purpose: 39* ============= 40*> 41*> \verbatim 42*> 43*> ZSTEMR computes selected eigenvalues and, optionally, eigenvectors 44*> of a real symmetric tridiagonal matrix T. Any such unreduced matrix has 45*> a well defined set of pairwise different real eigenvalues, the corresponding 46*> real eigenvectors are pairwise orthogonal. 47*> 48*> The spectrum may be computed either completely or partially by specifying 49*> either an interval (VL,VU] or a range of indices IL:IU for the desired 50*> eigenvalues. 51*> 52*> Depending on the number of desired eigenvalues, these are computed either 53*> by bisection or the dqds algorithm. Numerically orthogonal eigenvectors are 54*> computed by the use of various suitable L D L^T factorizations near clusters 55*> of close eigenvalues (referred to as RRRs, Relatively Robust 56*> Representations). An informal sketch of the algorithm follows. 57*> 58*> For each unreduced block (submatrix) of T, 59*> (a) Compute T - sigma I = L D L^T, so that L and D 60*> define all the wanted eigenvalues to high relative accuracy. 61*> This means that small relative changes in the entries of D and L 62*> cause only small relative changes in the eigenvalues and 63*> eigenvectors. The standard (unfactored) representation of the 64*> tridiagonal matrix T does not have this property in general. 65*> (b) Compute the eigenvalues to suitable accuracy. 66*> If the eigenvectors are desired, the algorithm attains full 67*> accuracy of the computed eigenvalues only right before 68*> the corresponding vectors have to be computed, see steps c) and d). 69*> (c) For each cluster of close eigenvalues, select a new 70*> shift close to the cluster, find a new factorization, and refine 71*> the shifted eigenvalues to suitable accuracy. 72*> (d) For each eigenvalue with a large enough relative separation compute 73*> the corresponding eigenvector by forming a rank revealing twisted 74*> factorization. Go back to (c) for any clusters that remain. 75*> 76*> For more details, see: 77*> - Inderjit S. Dhillon and Beresford N. Parlett: "Multiple representations 78*> to compute orthogonal eigenvectors of symmetric tridiagonal matrices," 79*> Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004. 80*> - Inderjit Dhillon and Beresford Parlett: "Orthogonal Eigenvectors and 81*> Relative Gaps," SIAM Journal on Matrix Analysis and Applications, Vol. 25, 82*> 2004. Also LAPACK Working Note 154. 83*> - Inderjit Dhillon: "A new O(n^2) algorithm for the symmetric 84*> tridiagonal eigenvalue/eigenvector problem", 85*> Computer Science Division Technical Report No. UCB/CSD-97-971, 86*> UC Berkeley, May 1997. 87*> 88*> Further Details 89*> 1.ZSTEMR works only on machines which follow IEEE-754 90*> floating-point standard in their handling of infinities and NaNs. 91*> This permits the use of efficient inner loops avoiding a check for 92*> zero divisors. 93*> 94*> 2. LAPACK routines can be used to reduce a complex Hermitean matrix to 95*> real symmetric tridiagonal form. 96*> 97*> (Any complex Hermitean tridiagonal matrix has real values on its diagonal 98*> and potentially complex numbers on its off-diagonals. By applying a 99*> similarity transform with an appropriate diagonal matrix 100*> diag(1,e^{i \phy_1}, ... , e^{i \phy_{n-1}}), the complex Hermitean 101*> matrix can be transformed into a real symmetric matrix and complex 102*> arithmetic can be entirely avoided.) 103*> 104*> While the eigenvectors of the real symmetric tridiagonal matrix are real, 105*> the eigenvectors of original complex Hermitean matrix have complex entries 106*> in general. 107*> Since LAPACK drivers overwrite the matrix data with the eigenvectors, 108*> ZSTEMR accepts complex workspace to facilitate interoperability 109*> with ZUNMTR or ZUPMTR. 110*> \endverbatim 111* 112* Arguments: 113* ========== 114* 115*> \param[in] JOBZ 116*> \verbatim 117*> JOBZ is CHARACTER*1 118*> = 'N': Compute eigenvalues only; 119*> = 'V': Compute eigenvalues and eigenvectors. 120*> \endverbatim 121*> 122*> \param[in] RANGE 123*> \verbatim 124*> RANGE is CHARACTER*1 125*> = 'A': all eigenvalues will be found. 126*> = 'V': all eigenvalues in the half-open interval (VL,VU] 127*> will be found. 128*> = 'I': the IL-th through IU-th eigenvalues will be found. 129*> \endverbatim 130*> 131*> \param[in] N 132*> \verbatim 133*> N is INTEGER 134*> The order of the matrix. N >= 0. 135*> \endverbatim 136*> 137*> \param[in,out] D 138*> \verbatim 139*> D is DOUBLE PRECISION array, dimension (N) 140*> On entry, the N diagonal elements of the tridiagonal matrix 141*> T. On exit, D is overwritten. 142*> \endverbatim 143*> 144*> \param[in,out] E 145*> \verbatim 146*> E is DOUBLE PRECISION array, dimension (N) 147*> On entry, the (N-1) subdiagonal elements of the tridiagonal 148*> matrix T in elements 1 to N-1 of E. E(N) need not be set on 149*> input, but is used internally as workspace. 150*> On exit, E is overwritten. 151*> \endverbatim 152*> 153*> \param[in] VL 154*> \verbatim 155*> VL is DOUBLE PRECISION 156*> \endverbatim 157*> 158*> \param[in] VU 159*> \verbatim 160*> VU is DOUBLE PRECISION 161*> 162*> If RANGE='V', the lower and upper bounds of the interval to 163*> be searched for eigenvalues. VL < VU. 164*> Not referenced if RANGE = 'A' or 'I'. 165*> \endverbatim 166*> 167*> \param[in] IL 168*> \verbatim 169*> IL is INTEGER 170*> \endverbatim 171*> 172*> \param[in] IU 173*> \verbatim 174*> IU is INTEGER 175*> 176*> If RANGE='I', the indices (in ascending order) of the 177*> smallest and largest eigenvalues to be returned. 178*> 1 <= IL <= IU <= N, if N > 0. 179*> Not referenced if RANGE = 'A' or 'V'. 180*> \endverbatim 181*> 182*> \param[out] M 183*> \verbatim 184*> M is INTEGER 185*> The total number of eigenvalues found. 0 <= M <= N. 186*> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1. 187*> \endverbatim 188*> 189*> \param[out] W 190*> \verbatim 191*> W is DOUBLE PRECISION array, dimension (N) 192*> The first M elements contain the selected eigenvalues in 193*> ascending order. 194*> \endverbatim 195*> 196*> \param[out] Z 197*> \verbatim 198*> Z is COMPLEX*16 array, dimension (LDZ, max(1,M) ) 199*> If JOBZ = 'V', and if INFO = 0, then the first M columns of Z 200*> contain the orthonormal eigenvectors of the matrix T 201*> corresponding to the selected eigenvalues, with the i-th 202*> column of Z holding the eigenvector associated with W(i). 203*> If JOBZ = 'N', then Z is not referenced. 204*> Note: the user must ensure that at least max(1,M) columns are 205*> supplied in the array Z; if RANGE = 'V', the exact value of M 206*> is not known in advance and can be computed with a workspace 207*> query by setting NZC = -1, see below. 208*> \endverbatim 209*> 210*> \param[in] LDZ 211*> \verbatim 212*> LDZ is INTEGER 213*> The leading dimension of the array Z. LDZ >= 1, and if 214*> JOBZ = 'V', then LDZ >= max(1,N). 215*> \endverbatim 216*> 217*> \param[in] NZC 218*> \verbatim 219*> NZC is INTEGER 220*> The number of eigenvectors to be held in the array Z. 221*> If RANGE = 'A', then NZC >= max(1,N). 222*> If RANGE = 'V', then NZC >= the number of eigenvalues in (VL,VU]. 223*> If RANGE = 'I', then NZC >= IU-IL+1. 224*> If NZC = -1, then a workspace query is assumed; the 225*> routine calculates the number of columns of the array Z that 226*> are needed to hold the eigenvectors. 227*> This value is returned as the first entry of the Z array, and 228*> no error message related to NZC is issued by XERBLA. 229*> \endverbatim 230*> 231*> \param[out] ISUPPZ 232*> \verbatim 233*> ISUPPZ is INTEGER ARRAY, dimension ( 2*max(1,M) ) 234*> The support of the eigenvectors in Z, i.e., the indices 235*> indicating the nonzero elements in Z. The i-th computed eigenvector 236*> is nonzero only in elements ISUPPZ( 2*i-1 ) through 237*> ISUPPZ( 2*i ). This is relevant in the case when the matrix 238*> is split. ISUPPZ is only accessed when JOBZ is 'V' and N > 0. 239*> \endverbatim 240*> 241*> \param[in,out] TRYRAC 242*> \verbatim 243*> TRYRAC is LOGICAL 244*> If TRYRAC.EQ..TRUE., indicates that the code should check whether 245*> the tridiagonal matrix defines its eigenvalues to high relative 246*> accuracy. If so, the code uses relative-accuracy preserving 247*> algorithms that might be (a bit) slower depending on the matrix. 248*> If the matrix does not define its eigenvalues to high relative 249*> accuracy, the code can uses possibly faster algorithms. 250*> If TRYRAC.EQ..FALSE., the code is not required to guarantee 251*> relatively accurate eigenvalues and can use the fastest possible 252*> techniques. 253*> On exit, a .TRUE. TRYRAC will be set to .FALSE. if the matrix 254*> does not define its eigenvalues to high relative accuracy. 255*> \endverbatim 256*> 257*> \param[out] WORK 258*> \verbatim 259*> WORK is DOUBLE PRECISION array, dimension (LWORK) 260*> On exit, if INFO = 0, WORK(1) returns the optimal 261*> (and minimal) LWORK. 262*> \endverbatim 263*> 264*> \param[in] LWORK 265*> \verbatim 266*> LWORK is INTEGER 267*> The dimension of the array WORK. LWORK >= max(1,18*N) 268*> if JOBZ = 'V', and LWORK >= max(1,12*N) if JOBZ = 'N'. 269*> If LWORK = -1, then a workspace query is assumed; the routine 270*> only calculates the optimal size of the WORK array, returns 271*> this value as the first entry of the WORK array, and no error 272*> message related to LWORK is issued by XERBLA. 273*> \endverbatim 274*> 275*> \param[out] IWORK 276*> \verbatim 277*> IWORK is INTEGER array, dimension (LIWORK) 278*> On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK. 279*> \endverbatim 280*> 281*> \param[in] LIWORK 282*> \verbatim 283*> LIWORK is INTEGER 284*> The dimension of the array IWORK. LIWORK >= max(1,10*N) 285*> if the eigenvectors are desired, and LIWORK >= max(1,8*N) 286*> if only the eigenvalues are to be computed. 287*> If LIWORK = -1, then a workspace query is assumed; the 288*> routine only calculates the optimal size of the IWORK array, 289*> returns this value as the first entry of the IWORK array, and 290*> no error message related to LIWORK is issued by XERBLA. 291*> \endverbatim 292*> 293*> \param[out] INFO 294*> \verbatim 295*> INFO is INTEGER 296*> On exit, INFO 297*> = 0: successful exit 298*> < 0: if INFO = -i, the i-th argument had an illegal value 299*> > 0: if INFO = 1X, internal error in DLARRE, 300*> if INFO = 2X, internal error in ZLARRV. 301*> Here, the digit X = ABS( IINFO ) < 10, where IINFO is 302*> the nonzero error code returned by DLARRE or 303*> ZLARRV, respectively. 304*> \endverbatim 305* 306* Authors: 307* ======== 308* 309*> \author Univ. of Tennessee 310*> \author Univ. of California Berkeley 311*> \author Univ. of Colorado Denver 312*> \author NAG Ltd. 313* 314*> \date November 2015 315* 316*> \ingroup complex16OTHERcomputational 317* 318*> \par Contributors: 319* ================== 320*> 321*> Beresford Parlett, University of California, Berkeley, USA \n 322*> Jim Demmel, University of California, Berkeley, USA \n 323*> Inderjit Dhillon, University of Texas, Austin, USA \n 324*> Osni Marques, LBNL/NERSC, USA \n 325*> Christof Voemel, University of California, Berkeley, USA \n 326* 327* ===================================================================== 328 SUBROUTINE ZSTEMR( JOBZ, RANGE, N, D, E, VL, VU, IL, IU, 329 $ M, W, Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, 330 $ IWORK, LIWORK, INFO ) 331* 332* -- LAPACK computational routine (version 3.6.0) -- 333* -- LAPACK is a software package provided by Univ. of Tennessee, -- 334* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 335* November 2015 336* 337* .. Scalar Arguments .. 338 CHARACTER JOBZ, RANGE 339 LOGICAL TRYRAC 340 INTEGER IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N 341 DOUBLE PRECISION VL, VU 342* .. 343* .. Array Arguments .. 344 INTEGER ISUPPZ( * ), IWORK( * ) 345 DOUBLE PRECISION D( * ), E( * ), W( * ), WORK( * ) 346 COMPLEX*16 Z( LDZ, * ) 347* .. 348* 349* ===================================================================== 350* 351* .. Parameters .. 352 DOUBLE PRECISION ZERO, ONE, FOUR, MINRGP 353 PARAMETER ( ZERO = 0.0D0, ONE = 1.0D0, 354 $ FOUR = 4.0D0, 355 $ MINRGP = 1.0D-3 ) 356* .. 357* .. Local Scalars .. 358 LOGICAL ALLEIG, INDEIG, LQUERY, VALEIG, WANTZ, ZQUERY 359 INTEGER I, IBEGIN, IEND, IFIRST, IIL, IINDBL, IINDW, 360 $ IINDWK, IINFO, IINSPL, IIU, ILAST, IN, INDD, 361 $ INDE2, INDERR, INDGP, INDGRS, INDWRK, ITMP, 362 $ ITMP2, J, JBLK, JJ, LIWMIN, LWMIN, NSPLIT, 363 $ NZCMIN, OFFSET, WBEGIN, WEND 364 DOUBLE PRECISION BIGNUM, CS, EPS, PIVMIN, R1, R2, RMAX, RMIN, 365 $ RTOL1, RTOL2, SAFMIN, SCALE, SMLNUM, SN, 366 $ THRESH, TMP, TNRM, WL, WU 367* .. 368* .. 369* .. External Functions .. 370 LOGICAL LSAME 371 DOUBLE PRECISION DLAMCH, DLANST 372 EXTERNAL LSAME, DLAMCH, DLANST 373* .. 374* .. External Subroutines .. 375 EXTERNAL DCOPY, DLAE2, DLAEV2, DLARRC, DLARRE, DLARRJ, 376 $ DLARRR, DLASRT, DSCAL, XERBLA, ZLARRV, ZSWAP 377* .. 378* .. Intrinsic Functions .. 379 INTRINSIC MAX, MIN, SQRT 380 381 382* .. 383* .. Executable Statements .. 384* 385* Test the input parameters. 386* 387 WANTZ = LSAME( JOBZ, 'V' ) 388 ALLEIG = LSAME( RANGE, 'A' ) 389 VALEIG = LSAME( RANGE, 'V' ) 390 INDEIG = LSAME( RANGE, 'I' ) 391* 392 LQUERY = ( ( LWORK.EQ.-1 ).OR.( LIWORK.EQ.-1 ) ) 393 ZQUERY = ( NZC.EQ.-1 ) 394 395* DSTEMR needs WORK of size 6*N, IWORK of size 3*N. 396* In addition, DLARRE needs WORK of size 6*N, IWORK of size 5*N. 397* Furthermore, ZLARRV needs WORK of size 12*N, IWORK of size 7*N. 398 IF( WANTZ ) THEN 399 LWMIN = 18*N 400 LIWMIN = 10*N 401 ELSE 402* need less workspace if only the eigenvalues are wanted 403 LWMIN = 12*N 404 LIWMIN = 8*N 405 ENDIF 406 407 WL = ZERO 408 WU = ZERO 409 IIL = 0 410 IIU = 0 411 NSPLIT = 0 412 413 IF( VALEIG ) THEN 414* We do not reference VL, VU in the cases RANGE = 'I','A' 415* The interval (WL, WU] contains all the wanted eigenvalues. 416* It is either given by the user or computed in DLARRE. 417 WL = VL 418 WU = VU 419 ELSEIF( INDEIG ) THEN 420* We do not reference IL, IU in the cases RANGE = 'V','A' 421 IIL = IL 422 IIU = IU 423 ENDIF 424* 425 INFO = 0 426 IF( .NOT.( WANTZ .OR. LSAME( JOBZ, 'N' ) ) ) THEN 427 INFO = -1 428 ELSE IF( .NOT.( ALLEIG .OR. VALEIG .OR. INDEIG ) ) THEN 429 INFO = -2 430 ELSE IF( N.LT.0 ) THEN 431 INFO = -3 432 ELSE IF( VALEIG .AND. N.GT.0 .AND. WU.LE.WL ) THEN 433 INFO = -7 434 ELSE IF( INDEIG .AND. ( IIL.LT.1 .OR. IIL.GT.N ) ) THEN 435 INFO = -8 436 ELSE IF( INDEIG .AND. ( IIU.LT.IIL .OR. IIU.GT.N ) ) THEN 437 INFO = -9 438 ELSE IF( LDZ.LT.1 .OR. ( WANTZ .AND. LDZ.LT.N ) ) THEN 439 INFO = -13 440 ELSE IF( LWORK.LT.LWMIN .AND. .NOT.LQUERY ) THEN 441 INFO = -17 442 ELSE IF( LIWORK.LT.LIWMIN .AND. .NOT.LQUERY ) THEN 443 INFO = -19 444 END IF 445* 446* Get machine constants. 447* 448 SAFMIN = DLAMCH( 'Safe minimum' ) 449 EPS = DLAMCH( 'Precision' ) 450 SMLNUM = SAFMIN / EPS 451 BIGNUM = ONE / SMLNUM 452 RMIN = SQRT( SMLNUM ) 453 RMAX = MIN( SQRT( BIGNUM ), ONE / SQRT( SQRT( SAFMIN ) ) ) 454* 455 IF( INFO.EQ.0 ) THEN 456 WORK( 1 ) = LWMIN 457 IWORK( 1 ) = LIWMIN 458* 459 IF( WANTZ .AND. ALLEIG ) THEN 460 NZCMIN = N 461 ELSE IF( WANTZ .AND. VALEIG ) THEN 462 CALL DLARRC( 'T', N, VL, VU, D, E, SAFMIN, 463 $ NZCMIN, ITMP, ITMP2, INFO ) 464 ELSE IF( WANTZ .AND. INDEIG ) THEN 465 NZCMIN = IIU-IIL+1 466 ELSE 467* WANTZ .EQ. FALSE. 468 NZCMIN = 0 469 ENDIF 470 IF( ZQUERY .AND. INFO.EQ.0 ) THEN 471 Z( 1,1 ) = NZCMIN 472 ELSE IF( NZC.LT.NZCMIN .AND. .NOT.ZQUERY ) THEN 473 INFO = -14 474 END IF 475 END IF 476 477 IF( INFO.NE.0 ) THEN 478* 479 CALL XERBLA( 'ZSTEMR', -INFO ) 480* 481 RETURN 482 ELSE IF( LQUERY .OR. ZQUERY ) THEN 483 RETURN 484 END IF 485* 486* Handle N = 0, 1, and 2 cases immediately 487* 488 M = 0 489 IF( N.EQ.0 ) 490 $ RETURN 491* 492 IF( N.EQ.1 ) THEN 493 IF( ALLEIG .OR. INDEIG ) THEN 494 M = 1 495 W( 1 ) = D( 1 ) 496 ELSE 497 IF( WL.LT.D( 1 ) .AND. WU.GE.D( 1 ) ) THEN 498 M = 1 499 W( 1 ) = D( 1 ) 500 END IF 501 END IF 502 IF( WANTZ.AND.(.NOT.ZQUERY) ) THEN 503 Z( 1, 1 ) = ONE 504 ISUPPZ(1) = 1 505 ISUPPZ(2) = 1 506 END IF 507 RETURN 508 END IF 509* 510 IF( N.EQ.2 ) THEN 511 IF( .NOT.WANTZ ) THEN 512 CALL DLAE2( D(1), E(1), D(2), R1, R2 ) 513 ELSE IF( WANTZ.AND.(.NOT.ZQUERY) ) THEN 514 CALL DLAEV2( D(1), E(1), D(2), R1, R2, CS, SN ) 515 END IF 516 IF( ALLEIG.OR. 517 $ (VALEIG.AND.(R2.GT.WL).AND. 518 $ (R2.LE.WU)).OR. 519 $ (INDEIG.AND.(IIL.EQ.1)) ) THEN 520 M = M+1 521 W( M ) = R2 522 IF( WANTZ.AND.(.NOT.ZQUERY) ) THEN 523 Z( 1, M ) = -SN 524 Z( 2, M ) = CS 525* Note: At most one of SN and CS can be zero. 526 IF (SN.NE.ZERO) THEN 527 IF (CS.NE.ZERO) THEN 528 ISUPPZ(2*M-1) = 1 529 ISUPPZ(2*M) = 2 530 ELSE 531 ISUPPZ(2*M-1) = 1 532 ISUPPZ(2*M) = 1 533 END IF 534 ELSE 535 ISUPPZ(2*M-1) = 2 536 ISUPPZ(2*M) = 2 537 END IF 538 ENDIF 539 ENDIF 540 IF( ALLEIG.OR. 541 $ (VALEIG.AND.(R1.GT.WL).AND. 542 $ (R1.LE.WU)).OR. 543 $ (INDEIG.AND.(IIU.EQ.2)) ) THEN 544 M = M+1 545 W( M ) = R1 546 IF( WANTZ.AND.(.NOT.ZQUERY) ) THEN 547 Z( 1, M ) = CS 548 Z( 2, M ) = SN 549* Note: At most one of SN and CS can be zero. 550 IF (SN.NE.ZERO) THEN 551 IF (CS.NE.ZERO) THEN 552 ISUPPZ(2*M-1) = 1 553 ISUPPZ(2*M) = 2 554 ELSE 555 ISUPPZ(2*M-1) = 1 556 ISUPPZ(2*M) = 1 557 END IF 558 ELSE 559 ISUPPZ(2*M-1) = 2 560 ISUPPZ(2*M) = 2 561 END IF 562 ENDIF 563 ENDIF 564 ELSE 565 566* Continue with general N 567 568 INDGRS = 1 569 INDERR = 2*N + 1 570 INDGP = 3*N + 1 571 INDD = 4*N + 1 572 INDE2 = 5*N + 1 573 INDWRK = 6*N + 1 574* 575 IINSPL = 1 576 IINDBL = N + 1 577 IINDW = 2*N + 1 578 IINDWK = 3*N + 1 579* 580* Scale matrix to allowable range, if necessary. 581* The allowable range is related to the PIVMIN parameter; see the 582* comments in DLARRD. The preference for scaling small values 583* up is heuristic; we expect users' matrices not to be close to the 584* RMAX threshold. 585* 586 SCALE = ONE 587 TNRM = DLANST( 'M', N, D, E ) 588 IF( TNRM.GT.ZERO .AND. TNRM.LT.RMIN ) THEN 589 SCALE = RMIN / TNRM 590 ELSE IF( TNRM.GT.RMAX ) THEN 591 SCALE = RMAX / TNRM 592 END IF 593 IF( SCALE.NE.ONE ) THEN 594 CALL DSCAL( N, SCALE, D, 1 ) 595 CALL DSCAL( N-1, SCALE, E, 1 ) 596 TNRM = TNRM*SCALE 597 IF( VALEIG ) THEN 598* If eigenvalues in interval have to be found, 599* scale (WL, WU] accordingly 600 WL = WL*SCALE 601 WU = WU*SCALE 602 ENDIF 603 END IF 604* 605* Compute the desired eigenvalues of the tridiagonal after splitting 606* into smaller subblocks if the corresponding off-diagonal elements 607* are small 608* THRESH is the splitting parameter for DLARRE 609* A negative THRESH forces the old splitting criterion based on the 610* size of the off-diagonal. A positive THRESH switches to splitting 611* which preserves relative accuracy. 612* 613 IF( TRYRAC ) THEN 614* Test whether the matrix warrants the more expensive relative approach. 615 CALL DLARRR( N, D, E, IINFO ) 616 ELSE 617* The user does not care about relative accurately eigenvalues 618 IINFO = -1 619 ENDIF 620* Set the splitting criterion 621 IF (IINFO.EQ.0) THEN 622 THRESH = EPS 623 ELSE 624 THRESH = -EPS 625* relative accuracy is desired but T does not guarantee it 626 TRYRAC = .FALSE. 627 ENDIF 628* 629 IF( TRYRAC ) THEN 630* Copy original diagonal, needed to guarantee relative accuracy 631 CALL DCOPY(N,D,1,WORK(INDD),1) 632 ENDIF 633* Store the squares of the offdiagonal values of T 634 DO 5 J = 1, N-1 635 WORK( INDE2+J-1 ) = E(J)**2 636 5 CONTINUE 637 638* Set the tolerance parameters for bisection 639 IF( .NOT.WANTZ ) THEN 640* DLARRE computes the eigenvalues to full precision. 641 RTOL1 = FOUR * EPS 642 RTOL2 = FOUR * EPS 643 ELSE 644* DLARRE computes the eigenvalues to less than full precision. 645* ZLARRV will refine the eigenvalue approximations, and we only 646* need less accurate initial bisection in DLARRE. 647* Note: these settings do only affect the subset case and DLARRE 648 RTOL1 = SQRT(EPS) 649 RTOL2 = MAX( SQRT(EPS)*5.0D-3, FOUR * EPS ) 650 ENDIF 651 CALL DLARRE( RANGE, N, WL, WU, IIL, IIU, D, E, 652 $ WORK(INDE2), RTOL1, RTOL2, THRESH, NSPLIT, 653 $ IWORK( IINSPL ), M, W, WORK( INDERR ), 654 $ WORK( INDGP ), IWORK( IINDBL ), 655 $ IWORK( IINDW ), WORK( INDGRS ), PIVMIN, 656 $ WORK( INDWRK ), IWORK( IINDWK ), IINFO ) 657 IF( IINFO.NE.0 ) THEN 658 INFO = 10 + ABS( IINFO ) 659 RETURN 660 END IF 661* Note that if RANGE .NE. 'V', DLARRE computes bounds on the desired 662* part of the spectrum. All desired eigenvalues are contained in 663* (WL,WU] 664 665 666 IF( WANTZ ) THEN 667* 668* Compute the desired eigenvectors corresponding to the computed 669* eigenvalues 670* 671 CALL ZLARRV( N, WL, WU, D, E, 672 $ PIVMIN, IWORK( IINSPL ), M, 673 $ 1, M, MINRGP, RTOL1, RTOL2, 674 $ W, WORK( INDERR ), WORK( INDGP ), IWORK( IINDBL ), 675 $ IWORK( IINDW ), WORK( INDGRS ), Z, LDZ, 676 $ ISUPPZ, WORK( INDWRK ), IWORK( IINDWK ), IINFO ) 677 IF( IINFO.NE.0 ) THEN 678 INFO = 20 + ABS( IINFO ) 679 RETURN 680 END IF 681 ELSE 682* DLARRE computes eigenvalues of the (shifted) root representation 683* ZLARRV returns the eigenvalues of the unshifted matrix. 684* However, if the eigenvectors are not desired by the user, we need 685* to apply the corresponding shifts from DLARRE to obtain the 686* eigenvalues of the original matrix. 687 DO 20 J = 1, M 688 ITMP = IWORK( IINDBL+J-1 ) 689 W( J ) = W( J ) + E( IWORK( IINSPL+ITMP-1 ) ) 690 20 CONTINUE 691 END IF 692* 693 694 IF ( TRYRAC ) THEN 695* Refine computed eigenvalues so that they are relatively accurate 696* with respect to the original matrix T. 697 IBEGIN = 1 698 WBEGIN = 1 699 DO 39 JBLK = 1, IWORK( IINDBL+M-1 ) 700 IEND = IWORK( IINSPL+JBLK-1 ) 701 IN = IEND - IBEGIN + 1 702 WEND = WBEGIN - 1 703* check if any eigenvalues have to be refined in this block 704 36 CONTINUE 705 IF( WEND.LT.M ) THEN 706 IF( IWORK( IINDBL+WEND ).EQ.JBLK ) THEN 707 WEND = WEND + 1 708 GO TO 36 709 END IF 710 END IF 711 IF( WEND.LT.WBEGIN ) THEN 712 IBEGIN = IEND + 1 713 GO TO 39 714 END IF 715 716 OFFSET = IWORK(IINDW+WBEGIN-1)-1 717 IFIRST = IWORK(IINDW+WBEGIN-1) 718 ILAST = IWORK(IINDW+WEND-1) 719 RTOL2 = FOUR * EPS 720 CALL DLARRJ( IN, 721 $ WORK(INDD+IBEGIN-1), WORK(INDE2+IBEGIN-1), 722 $ IFIRST, ILAST, RTOL2, OFFSET, W(WBEGIN), 723 $ WORK( INDERR+WBEGIN-1 ), 724 $ WORK( INDWRK ), IWORK( IINDWK ), PIVMIN, 725 $ TNRM, IINFO ) 726 IBEGIN = IEND + 1 727 WBEGIN = WEND + 1 728 39 CONTINUE 729 ENDIF 730* 731* If matrix was scaled, then rescale eigenvalues appropriately. 732* 733 IF( SCALE.NE.ONE ) THEN 734 CALL DSCAL( M, ONE / SCALE, W, 1 ) 735 END IF 736 END IF 737* 738* If eigenvalues are not in increasing order, then sort them, 739* possibly along with eigenvectors. 740* 741 IF( NSPLIT.GT.1 .OR. N.EQ.2 ) THEN 742 IF( .NOT. WANTZ ) THEN 743 CALL DLASRT( 'I', M, W, IINFO ) 744 IF( IINFO.NE.0 ) THEN 745 INFO = 3 746 RETURN 747 END IF 748 ELSE 749 DO 60 J = 1, M - 1 750 I = 0 751 TMP = W( J ) 752 DO 50 JJ = J + 1, M 753 IF( W( JJ ).LT.TMP ) THEN 754 I = JJ 755 TMP = W( JJ ) 756 END IF 757 50 CONTINUE 758 IF( I.NE.0 ) THEN 759 W( I ) = W( J ) 760 W( J ) = TMP 761 IF( WANTZ ) THEN 762 CALL ZSWAP( N, Z( 1, I ), 1, Z( 1, J ), 1 ) 763 ITMP = ISUPPZ( 2*I-1 ) 764 ISUPPZ( 2*I-1 ) = ISUPPZ( 2*J-1 ) 765 ISUPPZ( 2*J-1 ) = ITMP 766 ITMP = ISUPPZ( 2*I ) 767 ISUPPZ( 2*I ) = ISUPPZ( 2*J ) 768 ISUPPZ( 2*J ) = ITMP 769 END IF 770 END IF 771 60 CONTINUE 772 END IF 773 ENDIF 774* 775* 776 WORK( 1 ) = LWMIN 777 IWORK( 1 ) = LIWMIN 778 RETURN 779* 780* End of ZSTEMR 781* 782 END 783