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

◆ cstemr()

subroutine cstemr ( character jobz,
character range,
integer n,
real, dimension( * ) d,
real, dimension( * ) e,
real vl,
real vu,
integer il,
integer iu,
integer m,
real, dimension( * ) w,
complex, dimension( ldz, * ) z,
integer ldz,
integer nzc,
integer, dimension( * ) isuppz,
logical tryrac,
real, dimension( * ) work,
integer lwork,
integer, dimension( * ) iwork,
integer liwork,
integer info )

CSTEMR

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

Purpose:
!>
!> CSTEMR computes selected eigenvalues and, optionally, eigenvectors
!> of a real symmetric tridiagonal matrix T. Any such unreduced matrix has
!> a well defined set of pairwise different real eigenvalues, the corresponding
!> real eigenvectors are pairwise orthogonal.
!>
!> The spectrum may be computed either completely or partially by specifying
!> either an interval (VL,VU] or a range of indices IL:IU for the desired
!> eigenvalues.
!>
!> Depending on the number of desired eigenvalues, these are computed either
!> by bisection or the dqds algorithm. Numerically orthogonal eigenvectors are
!> computed by the use of various suitable L D L^T factorizations near clusters
!> of close eigenvalues (referred to as RRRs, Relatively Robust
!> Representations). An informal sketch of the algorithm follows.
!>
!> For each unreduced block (submatrix) of T,
!>    (a) Compute T - sigma I  = L D L^T, so that L and D
!>        define all the wanted eigenvalues to high relative accuracy.
!>        This means that small relative changes in the entries of D and L
!>        cause only small relative changes in the eigenvalues and
!>        eigenvectors. The standard (unfactored) representation of the
!>        tridiagonal matrix T does not have this property in general.
!>    (b) Compute the eigenvalues to suitable accuracy.
!>        If the eigenvectors are desired, the algorithm attains full
!>        accuracy of the computed eigenvalues only right before
!>        the corresponding vectors have to be computed, see steps c) and d).
!>    (c) For each cluster of close eigenvalues, select a new
!>        shift close to the cluster, find a new factorization, and refine
!>        the shifted eigenvalues to suitable accuracy.
!>    (d) For each eigenvalue with a large enough relative separation compute
!>        the corresponding eigenvector by forming a rank revealing twisted
!>        factorization. Go back to (c) for any clusters that remain.
!>
!> For more details, see:
!> - Inderjit S. Dhillon and Beresford N. Parlett: 
!>   Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004.
!> - Inderjit Dhillon and Beresford Parlett:  SIAM Journal on Matrix Analysis and Applications, Vol. 25,
!>   2004.  Also LAPACK Working Note 154.
!> - Inderjit Dhillon: ,
!>   Computer Science Division Technical Report No. UCB/CSD-97-971,
!>   UC Berkeley, May 1997.
!>
!> Further Details
!> 1.CSTEMR works only on machines which follow IEEE-754
!> floating-point standard in their handling of infinities and NaNs.
!> This permits the use of efficient inner loops avoiding a check for
!> zero divisors.
!>
!> 2. LAPACK routines can be used to reduce a complex Hermitean matrix to
!> real symmetric tridiagonal form.
!>
!> (Any complex Hermitean tridiagonal matrix has real values on its diagonal
!> and potentially complex numbers on its off-diagonals. By applying a
!> similarity transform with an appropriate diagonal matrix
!> diag(1,e^{i \phy_1}, ... , e^{i \phy_{n-1}}), the complex Hermitean
!> matrix can be transformed into a real symmetric matrix and complex
!> arithmetic can be entirely avoided.)
!>
!> While the eigenvectors of the real symmetric tridiagonal matrix are real,
!> the eigenvectors of original complex Hermitean matrix have complex entries
!> in general.
!> Since LAPACK drivers overwrite the matrix data with the eigenvectors,
!> CSTEMR accepts complex workspace to facilitate interoperability
!> with CUNMTR or CUPMTR.
!> 
Parameters
[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]N
!>          N is INTEGER
!>          The order of the matrix.  N >= 0.
!> 
[in,out]D
!>          D is REAL array, dimension (N)
!>          On entry, the N diagonal elements of the tridiagonal matrix
!>          T. On exit, D is overwritten.
!> 
[in,out]E
!>          E is REAL array, dimension (N)
!>          On entry, the (N-1) subdiagonal elements of the tridiagonal
!>          matrix T in elements 1 to N-1 of E. E(N) need not be set on
!>          input, but is used internally as workspace.
!>          On exit, E is overwritten.
!> 
[in]VL
!>          VL is REAL
!>
!>          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 REAL
!>
!>          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.
!>          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.
!>          Not referenced if RANGE = 'A' or 'V'.
!> 
[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 REAL array, dimension (N)
!>          The first M elements contain the selected eigenvalues in
!>          ascending order.
!> 
[out]Z
!>          Z is COMPLEX array, dimension (LDZ, max(1,M) )
!>          If JOBZ = 'V', and if INFO = 0, then the first M columns of Z
!>          contain the orthonormal eigenvectors of the matrix T
!>          corresponding to the selected eigenvalues, with the i-th
!>          column of Z holding the eigenvector associated with W(i).
!>          If JOBZ = 'N', then Z is not referenced.
!>          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 can be computed with a workspace
!>          query by setting NZC = -1, see below.
!> 
[in]LDZ
!>          LDZ is INTEGER
!>          The leading dimension of the array Z.  LDZ >= 1, and if
!>          JOBZ = 'V', then LDZ >= max(1,N).
!> 
[in]NZC
!>          NZC is INTEGER
!>          The number of eigenvectors to be held in the array Z.
!>          If RANGE = 'A', then NZC >= max(1,N).
!>          If RANGE = 'V', then NZC >= the number of eigenvalues in (VL,VU].
!>          If RANGE = 'I', then NZC >= IU-IL+1.
!>          If NZC = -1, then a workspace query is assumed; the
!>          routine calculates the number of columns of the array Z that
!>          are needed to hold the eigenvectors.
!>          This value is returned as the first entry of the Z array, and
!>          no error message related to NZC is issued by XERBLA.
!> 
[out]ISUPPZ
!>          ISUPPZ is INTEGER array, dimension ( 2*max(1,M) )
!>          The support of the eigenvectors in Z, i.e., the indices
!>          indicating the nonzero elements in Z. The i-th computed eigenvector
!>          is nonzero only in elements ISUPPZ( 2*i-1 ) through
!>          ISUPPZ( 2*i ). This is relevant in the case when the matrix
!>          is split. ISUPPZ is only accessed when JOBZ is 'V' and N > 0.
!> 
[in,out]TRYRAC
!>          TRYRAC is LOGICAL
!>          If TRYRAC = .TRUE., indicates that the code should check whether
!>          the tridiagonal matrix defines its eigenvalues to high relative
!>          accuracy.  If so, the code uses relative-accuracy preserving
!>          algorithms that might be (a bit) slower depending on the matrix.
!>          If the matrix does not define its eigenvalues to high relative
!>          accuracy, the code can uses possibly faster algorithms.
!>          If TRYRAC = .FALSE., the code is not required to guarantee
!>          relatively accurate eigenvalues and can use the fastest possible
!>          techniques.
!>          On exit, a .TRUE. TRYRAC will be set to .FALSE. if the matrix
!>          does not define its eigenvalues to high relative accuracy.
!> 
[out]WORK
!>          WORK is REAL array, dimension (LWORK)
!>          On exit, if INFO = 0, WORK(1) returns the optimal
!>          (and minimal) LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The dimension of the array WORK. LWORK >= max(1,18*N)
!>          if JOBZ = 'V', and LWORK >= max(1,12*N) if JOBZ = 'N'.
!>          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]IWORK
!>          IWORK is INTEGER array, dimension (LIWORK)
!>          On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK.
!> 
[in]LIWORK
!>          LIWORK is INTEGER
!>          The dimension of the array IWORK.  LIWORK >= max(1,10*N)
!>          if the eigenvectors are desired, and LIWORK >= max(1,8*N)
!>          if only the eigenvalues are to be computed.
!>          If LIWORK = -1, then a workspace query is assumed; the
!>          routine only calculates the optimal size of the IWORK array,
!>          returns this value as the first entry of the IWORK array, and
!>          no error message related to LIWORK is issued by XERBLA.
!> 
[out]INFO
!>          INFO is INTEGER
!>          On exit, INFO
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value
!>          > 0:  if INFO = 1X, internal error in SLARRE,
!>                if INFO = 2X, internal error in CLARRV.
!>                Here, the digit X = ABS( IINFO ) < 10, where IINFO is
!>                the nonzero error code returned by SLARRE or
!>                CLARRV, respectively.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Beresford Parlett, University of California, Berkeley, USA
Jim Demmel, University of California, Berkeley, USA
Inderjit Dhillon, University of Texas, Austin, USA
Osni Marques, LBNL/NERSC, USA
Christof Voemel, University of California, Berkeley, USA
Aravindh Krishnamoorthy, FAU, Erlangen, Germany

Definition at line 334 of file cstemr.f.

337*
338* -- LAPACK computational routine --
339* -- LAPACK is a software package provided by Univ. of Tennessee, --
340* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
341*
342* .. Scalar Arguments ..
343 CHARACTER JOBZ, RANGE
344 LOGICAL TRYRAC
345 INTEGER IL, INFO, IU, LDZ, NZC, LIWORK, LWORK, M, N
346 REAL VL, VU
347* ..
348* .. Array Arguments ..
349 INTEGER ISUPPZ( * ), IWORK( * )
350 REAL D( * ), E( * ), W( * ), WORK( * )
351 COMPLEX Z( LDZ, * )
352* ..
353*
354* =====================================================================
355*
356* .. Parameters ..
357 REAL ZERO, ONE, FOUR, MINRGP
358 parameter( zero = 0.0e0, one = 1.0e0,
359 $ four = 4.0e0,
360 $ minrgp = 3.0e-3 )
361* ..
362* .. Local Scalars ..
363 LOGICAL ALLEIG, INDEIG, LQUERY, VALEIG, WANTZ, ZQUERY,
364 $ LAESWAP
365 INTEGER I, IBEGIN, IEND, IFIRST, IIL, IINDBL, IINDW,
366 $ IINDWK, IINFO, IINSPL, IIU, ILAST, IN, INDD,
367 $ INDE2, INDERR, INDGP, INDGRS, INDWRK, ITMP,
368 $ ITMP2, J, JBLK, JJ, LIWMIN, LWMIN, NSPLIT,
369 $ NZCMIN, OFFSET, WBEGIN, WEND
370 REAL BIGNUM, CS, EPS, PIVMIN, R1, R2, RMAX, RMIN,
371 $ RTOL1, RTOL2, SAFMIN, SCALE, SMLNUM, SN,
372 $ THRESH, TMP, TNRM, WL, WU
373* ..
374* ..
375* .. External Functions ..
376 LOGICAL LSAME
377 REAL SLAMCH, SLANST, SROUNDUP_LWORK
378 EXTERNAL lsame, slamch, slanst,
380* ..
381* .. External Subroutines ..
382 EXTERNAL clarrv, cswap, scopy, slae2, slaev2,
383 $ slarrc,
385* ..
386* .. Intrinsic Functions ..
387 INTRINSIC max, min, sqrt
388
389
390* ..
391* .. Executable Statements ..
392*
393* Test the input parameters.
394*
395 wantz = lsame( jobz, 'V' )
396 alleig = lsame( range, 'A' )
397 valeig = lsame( range, 'V' )
398 indeig = lsame( range, 'I' )
399*
400 lquery = ( ( lwork.EQ.-1 ).OR.( liwork.EQ.-1 ) )
401 zquery = ( nzc.EQ.-1 )
402 laeswap = .false.
403
404* SSTEMR needs WORK of size 6*N, IWORK of size 3*N.
405* In addition, SLARRE needs WORK of size 6*N, IWORK of size 5*N.
406* Furthermore, CLARRV needs WORK of size 12*N, IWORK of size 7*N.
407 IF( wantz ) THEN
408 lwmin = 18*n
409 liwmin = 10*n
410 ELSE
411* need less workspace if only the eigenvalues are wanted
412 lwmin = 12*n
413 liwmin = 8*n
414 ENDIF
415
416 wl = zero
417 wu = zero
418 iil = 0
419 iiu = 0
420 nsplit = 0
421
422 IF( valeig ) THEN
423* We do not reference VL, VU in the cases RANGE = 'I','A'
424* The interval (WL, WU] contains all the wanted eigenvalues.
425* It is either given by the user or computed in SLARRE.
426 wl = vl
427 wu = vu
428 ELSEIF( indeig ) THEN
429* We do not reference IL, IU in the cases RANGE = 'V','A'
430 iil = il
431 iiu = iu
432 ENDIF
433*
434 info = 0
435 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
436 info = -1
437 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
438 info = -2
439 ELSE IF( n.LT.0 ) THEN
440 info = -3
441 ELSE IF( valeig .AND. n.GT.0 .AND. wu.LE.wl ) THEN
442 info = -7
443 ELSE IF( indeig .AND. ( iil.LT.1 .OR. iil.GT.n ) ) THEN
444 info = -8
445 ELSE IF( indeig .AND. ( iiu.LT.iil .OR. iiu.GT.n ) ) THEN
446 info = -9
447 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
448 info = -13
449 ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
450 info = -17
451 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
452 info = -19
453 END IF
454*
455* Get machine constants.
456*
457 safmin = slamch( 'Safe minimum' )
458 eps = slamch( 'Precision' )
459 smlnum = safmin / eps
460 bignum = one / smlnum
461 rmin = sqrt( smlnum )
462 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
463*
464 IF( info.EQ.0 ) THEN
465 work( 1 ) = sroundup_lwork(lwmin)
466 iwork( 1 ) = liwmin
467*
468 IF( wantz .AND. alleig ) THEN
469 nzcmin = n
470 ELSE IF( wantz .AND. valeig ) THEN
471 CALL slarrc( 'T', n, vl, vu, d, e, safmin,
472 $ nzcmin, itmp, itmp2, info )
473 ELSE IF( wantz .AND. indeig ) THEN
474 nzcmin = iiu-iil+1
475 ELSE
476* WANTZ .EQ. FALSE.
477 nzcmin = 0
478 ENDIF
479 IF( zquery .AND. info.EQ.0 ) THEN
480 z( 1,1 ) = cmplx( nzcmin )
481 ELSE IF( nzc.LT.nzcmin .AND. .NOT.zquery ) THEN
482 info = -14
483 END IF
484 END IF
485
486 IF( info.NE.0 ) THEN
487*
488 CALL xerbla( 'CSTEMR', -info )
489*
490 RETURN
491 ELSE IF( lquery .OR. zquery ) THEN
492 RETURN
493 END IF
494*
495* Handle N = 0, 1, and 2 cases immediately
496*
497 m = 0
498 IF( n.EQ.0 )
499 $ RETURN
500*
501 IF( n.EQ.1 ) THEN
502 IF( alleig .OR. indeig ) THEN
503 m = 1
504 w( 1 ) = d( 1 )
505 ELSE
506 IF( wl.LT.d( 1 ) .AND. wu.GE.d( 1 ) ) THEN
507 m = 1
508 w( 1 ) = d( 1 )
509 END IF
510 END IF
511 IF( wantz.AND.(.NOT.zquery) ) THEN
512 z( 1, 1 ) = one
513 isuppz(1) = 1
514 isuppz(2) = 1
515 END IF
516 RETURN
517 END IF
518*
519 IF( n.EQ.2 ) THEN
520 IF( .NOT.wantz ) THEN
521 CALL slae2( d(1), e(1), d(2), r1, r2 )
522 ELSE IF( wantz.AND.(.NOT.zquery) ) THEN
523 CALL slaev2( d(1), e(1), d(2), r1, r2, cs, sn )
524 END IF
525* D/S/LAE2 and D/S/LAEV2 outputs satisfy |R1| >= |R2|. However,
526* the following code requires R1 >= R2. Hence, we correct
527* the order of R1, R2, CS, SN if R1 < R2 before further processing.
528 IF( r1.LT.r2 ) THEN
529 e(2) = r1
530 r1 = r2
531 r2 = e(2)
532 laeswap = .true.
533 ENDIF
534 IF( alleig.OR.
535 $ (valeig.AND.(r2.GT.wl).AND.
536 $ (r2.LE.wu)).OR.
537 $ (indeig.AND.(iil.EQ.1)) ) THEN
538 m = m+1
539 w( m ) = r2
540 IF( wantz.AND.(.NOT.zquery) ) THEN
541 IF( laeswap ) THEN
542 z( 1, m ) = cs
543 z( 2, m ) = sn
544 ELSE
545 z( 1, m ) = -sn
546 z( 2, m ) = cs
547 ENDIF
548* Note: At most one of SN and CS can be zero.
549 IF (sn.NE.zero) THEN
550 IF (cs.NE.zero) THEN
551 isuppz(2*m-1) = 1
552 isuppz(2*m) = 2
553 ELSE
554 isuppz(2*m-1) = 1
555 isuppz(2*m) = 1
556 END IF
557 ELSE
558 isuppz(2*m-1) = 2
559 isuppz(2*m) = 2
560 END IF
561 ENDIF
562 ENDIF
563 IF( alleig.OR.
564 $ (valeig.AND.(r1.GT.wl).AND.
565 $ (r1.LE.wu)).OR.
566 $ (indeig.AND.(iiu.EQ.2)) ) THEN
567 m = m+1
568 w( m ) = r1
569 IF( wantz.AND.(.NOT.zquery) ) THEN
570 IF( laeswap ) THEN
571 z( 1, m ) = -sn
572 z( 2, m ) = cs
573 ELSE
574 z( 1, m ) = cs
575 z( 2, m ) = sn
576 ENDIF
577* Note: At most one of SN and CS can be zero.
578 IF (sn.NE.zero) THEN
579 IF (cs.NE.zero) THEN
580 isuppz(2*m-1) = 1
581 isuppz(2*m) = 2
582 ELSE
583 isuppz(2*m-1) = 1
584 isuppz(2*m) = 1
585 END IF
586 ELSE
587 isuppz(2*m-1) = 2
588 isuppz(2*m) = 2
589 END IF
590 ENDIF
591 ENDIF
592 ELSE
593
594* Continue with general N
595
596 indgrs = 1
597 inderr = 2*n + 1
598 indgp = 3*n + 1
599 indd = 4*n + 1
600 inde2 = 5*n + 1
601 indwrk = 6*n + 1
602*
603 iinspl = 1
604 iindbl = n + 1
605 iindw = 2*n + 1
606 iindwk = 3*n + 1
607*
608* Scale matrix to allowable range, if necessary.
609* The allowable range is related to the PIVMIN parameter; see the
610* comments in SLARRD. The preference for scaling small values
611* up is heuristic; we expect users' matrices not to be close to the
612* RMAX threshold.
613*
614 scale = one
615 tnrm = slanst( 'M', n, d, e )
616 IF( tnrm.GT.zero .AND. tnrm.LT.rmin ) THEN
617 scale = rmin / tnrm
618 ELSE IF( tnrm.GT.rmax ) THEN
619 scale = rmax / tnrm
620 END IF
621 IF( scale.NE.one ) THEN
622 CALL sscal( n, scale, d, 1 )
623 CALL sscal( n-1, scale, e, 1 )
624 tnrm = tnrm*scale
625 IF( valeig ) THEN
626* If eigenvalues in interval have to be found,
627* scale (WL, WU] accordingly
628 wl = wl*scale
629 wu = wu*scale
630 ENDIF
631 END IF
632*
633* Compute the desired eigenvalues of the tridiagonal after splitting
634* into smaller subblocks if the corresponding off-diagonal elements
635* are small
636* THRESH is the splitting parameter for SLARRE
637* A negative THRESH forces the old splitting criterion based on the
638* size of the off-diagonal. A positive THRESH switches to splitting
639* which preserves relative accuracy.
640*
641 IF( tryrac ) THEN
642* Test whether the matrix warrants the more expensive relative approach.
643 CALL slarrr( n, d, e, iinfo )
644 ELSE
645* The user does not care about relative accurately eigenvalues
646 iinfo = -1
647 ENDIF
648* Set the splitting criterion
649 IF (iinfo.EQ.0) THEN
650 thresh = eps
651 ELSE
652 thresh = -eps
653* relative accuracy is desired but T does not guarantee it
654 tryrac = .false.
655 ENDIF
656*
657 IF( tryrac ) THEN
658* Copy original diagonal, needed to guarantee relative accuracy
659 CALL scopy(n,d,1,work(indd),1)
660 ENDIF
661* Store the squares of the offdiagonal values of T
662 DO 5 j = 1, n-1
663 work( inde2+j-1 ) = e(j)**2
664 5 CONTINUE
665
666* Set the tolerance parameters for bisection
667 IF( .NOT.wantz ) THEN
668* SLARRE computes the eigenvalues to full precision.
669 rtol1 = four * eps
670 rtol2 = four * eps
671 ELSE
672* SLARRE computes the eigenvalues to less than full precision.
673* CLARRV will refine the eigenvalue approximations, and we only
674* need less accurate initial bisection in SLARRE.
675* Note: these settings do only affect the subset case and SLARRE
676 rtol1 = max( sqrt(eps)*5.0e-2, four * eps )
677 rtol2 = max( sqrt(eps)*5.0e-3, four * eps )
678 ENDIF
679 CALL slarre( range, n, wl, wu, iil, iiu, d, e,
680 $ work(inde2), rtol1, rtol2, thresh, nsplit,
681 $ iwork( iinspl ), m, w, work( inderr ),
682 $ work( indgp ), iwork( iindbl ),
683 $ iwork( iindw ), work( indgrs ), pivmin,
684 $ work( indwrk ), iwork( iindwk ), iinfo )
685 IF( iinfo.NE.0 ) THEN
686 info = 10 + abs( iinfo )
687 RETURN
688 END IF
689* Note that if RANGE .NE. 'V', SLARRE computes bounds on the desired
690* part of the spectrum. All desired eigenvalues are contained in
691* (WL,WU]
692
693
694 IF( wantz ) THEN
695*
696* Compute the desired eigenvectors corresponding to the computed
697* eigenvalues
698*
699 CALL clarrv( n, wl, wu, d, e,
700 $ pivmin, iwork( iinspl ), m,
701 $ 1, m, minrgp, rtol1, rtol2,
702 $ w, work( inderr ), work( indgp ), iwork( iindbl ),
703 $ iwork( iindw ), work( indgrs ), z, ldz,
704 $ isuppz, work( indwrk ), iwork( iindwk ), iinfo )
705 IF( iinfo.NE.0 ) THEN
706 info = 20 + abs( iinfo )
707 RETURN
708 END IF
709 ELSE
710* SLARRE computes eigenvalues of the (shifted) root representation
711* CLARRV returns the eigenvalues of the unshifted matrix.
712* However, if the eigenvectors are not desired by the user, we need
713* to apply the corresponding shifts from SLARRE to obtain the
714* eigenvalues of the original matrix.
715 DO 20 j = 1, m
716 itmp = iwork( iindbl+j-1 )
717 w( j ) = w( j ) + e( iwork( iinspl+itmp-1 ) )
718 20 CONTINUE
719 END IF
720*
721
722 IF ( tryrac ) THEN
723* Refine computed eigenvalues so that they are relatively accurate
724* with respect to the original matrix T.
725 ibegin = 1
726 wbegin = 1
727 DO 39 jblk = 1, iwork( iindbl+m-1 )
728 iend = iwork( iinspl+jblk-1 )
729 in = iend - ibegin + 1
730 wend = wbegin - 1
731* check if any eigenvalues have to be refined in this block
732 36 CONTINUE
733 IF( wend.LT.m ) THEN
734 IF( iwork( iindbl+wend ).EQ.jblk ) THEN
735 wend = wend + 1
736 GO TO 36
737 END IF
738 END IF
739 IF( wend.LT.wbegin ) THEN
740 ibegin = iend + 1
741 GO TO 39
742 END IF
743
744 offset = iwork(iindw+wbegin-1)-1
745 ifirst = iwork(iindw+wbegin-1)
746 ilast = iwork(iindw+wend-1)
747 rtol2 = four * eps
748 CALL slarrj( in,
749 $ work(indd+ibegin-1), work(inde2+ibegin-1),
750 $ ifirst, ilast, rtol2, offset, w(wbegin),
751 $ work( inderr+wbegin-1 ),
752 $ work( indwrk ), iwork( iindwk ), pivmin,
753 $ tnrm, iinfo )
754 ibegin = iend + 1
755 wbegin = wend + 1
756 39 CONTINUE
757 ENDIF
758*
759* If matrix was scaled, then rescale eigenvalues appropriately.
760*
761 IF( scale.NE.one ) THEN
762 CALL sscal( m, one / scale, w, 1 )
763 END IF
764 END IF
765*
766* If eigenvalues are not in increasing order, then sort them,
767* possibly along with eigenvectors.
768*
769 IF( nsplit.GT.1 .OR. n.EQ.2 ) THEN
770 IF( .NOT. wantz ) THEN
771 CALL slasrt( 'I', m, w, iinfo )
772 IF( iinfo.NE.0 ) THEN
773 info = 3
774 RETURN
775 END IF
776 ELSE
777 DO 60 j = 1, m - 1
778 i = 0
779 tmp = w( j )
780 DO 50 jj = j + 1, m
781 IF( w( jj ).LT.tmp ) THEN
782 i = jj
783 tmp = w( jj )
784 END IF
785 50 CONTINUE
786 IF( i.NE.0 ) THEN
787 w( i ) = w( j )
788 w( j ) = tmp
789 IF( wantz ) THEN
790 CALL cswap( n, z( 1, i ), 1, z( 1, j ), 1 )
791 itmp = isuppz( 2*i-1 )
792 isuppz( 2*i-1 ) = isuppz( 2*j-1 )
793 isuppz( 2*j-1 ) = itmp
794 itmp = isuppz( 2*i )
795 isuppz( 2*i ) = isuppz( 2*j )
796 isuppz( 2*j ) = itmp
797 END IF
798 END IF
799 60 CONTINUE
800 END IF
801 ENDIF
802*
803*
804 work( 1 ) = sroundup_lwork(lwmin)
805 iwork( 1 ) = liwmin
806 RETURN
807*
808* End of CSTEMR
809*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine slae2(a, b, c, rt1, rt2)
SLAE2 computes the eigenvalues of a 2-by-2 symmetric matrix.
Definition slae2.f:100
subroutine slaev2(a, b, c, rt1, rt2, cs1, sn1)
SLAEV2 computes the eigenvalues and eigenvectors of a 2-by-2 symmetric/Hermitian matrix.
Definition slaev2.f:118
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function slanst(norm, n, d, e)
SLANST returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition slanst.f:98
subroutine slarrc(jobt, n, vl, vu, d, e, pivmin, eigcnt, lcnt, rcnt, info)
SLARRC computes the number of eigenvalues of the symmetric tridiagonal matrix.
Definition slarrc.f:135
subroutine slarre(range, n, vl, vu, il, iu, d, e, e2, rtol1, rtol2, spltol, nsplit, isplit, m, w, werr, wgap, iblock, indexw, gers, pivmin, work, iwork, info)
SLARRE given the tridiagonal matrix T, sets small off-diagonal elements to zero and for each unreduce...
Definition slarre.f:303
subroutine slarrj(n, d, e2, ifirst, ilast, rtol, offset, w, werr, work, iwork, pivmin, spdiam, info)
SLARRJ performs refinement of the initial estimates of the eigenvalues of the matrix T.
Definition slarrj.f:166
subroutine slarrr(n, d, e, info)
SLARRR performs tests to decide whether the symmetric tridiagonal matrix T warrants expensive computa...
Definition slarrr.f:92
subroutine clarrv(n, vl, vu, d, l, pivmin, isplit, m, dol, dou, minrgp, rtol1, rtol2, w, werr, wgap, iblock, indexw, gers, z, ldz, isuppz, work, iwork, info)
CLARRV computes the eigenvectors of the tridiagonal matrix T = L D LT given L, D and the eigenvalues ...
Definition clarrv.f:284
subroutine slasrt(id, n, d, info)
SLASRT sorts numbers in increasing or decreasing order.
Definition slasrt.f:86
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine cswap(n, cx, incx, cy, incy)
CSWAP
Definition cswap.f:81
Here is the call graph for this function:
Here is the caller graph for this function: