LAPACK  3.10.0
LAPACK: Linear Algebra PACKage

◆ cheevr()

subroutine cheevr ( character  JOBZ,
character  RANGE,
character  UPLO,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
real  VL,
real  VU,
integer  IL,
integer  IU,
real  ABSTOL,
integer  M,
real, dimension( * )  W,
complex, dimension( ldz, * )  Z,
integer  LDZ,
integer, dimension( * )  ISUPPZ,
complex, dimension( * )  WORK,
integer  LWORK,
real, dimension( * )  RWORK,
integer  LRWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

CHEEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices

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

Purpose:
 CHEEVR computes selected eigenvalues and, optionally, eigenvectors
 of a complex Hermitian matrix A.  Eigenvalues and eigenvectors can
 be selected by specifying either a range of values or a range of
 indices for the desired eigenvalues.

 CHEEVR first reduces the matrix A to tridiagonal form T with a call
 to CHETRD.  Then, whenever possible, CHEEVR calls CSTEMR to compute
 the eigenspectrum using Relatively Robust Representations.  CSTEMR
 computes eigenvalues by the dqds algorithm, while orthogonal
 eigenvectors are computed from various "good" L D L^T representations
 (also known as Relatively Robust Representations). Gram-Schmidt
 orthogonalization is avoided as far as possible. More specifically,
 the various steps of the algorithm are as 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.

 The desired accuracy of the output can be specified by the input
 parameter ABSTOL.

 For more details, see CSTEMR's documentation and:
 - Inderjit S. Dhillon and Beresford N. Parlett: "Multiple representations
   to compute orthogonal eigenvectors of symmetric tridiagonal matrices,"
   Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004.
 - Inderjit Dhillon and Beresford Parlett: "Orthogonal Eigenvectors and
   Relative Gaps," SIAM Journal on Matrix Analysis and Applications, Vol. 25,
   2004.  Also LAPACK Working Note 154.
 - Inderjit Dhillon: "A new O(n^2) algorithm for the symmetric
   tridiagonal eigenvalue/eigenvector problem",
   Computer Science Division Technical Report No. UCB/CSD-97-971,
   UC Berkeley, May 1997.


 Note 1 : CHEEVR calls CSTEMR when the full spectrum is requested
 on machines which conform to the ieee-754 floating point standard.
 CHEEVR calls SSTEBZ and CSTEIN on non-ieee machines and
 when partial spectrum requests are made.

 Normal execution of CSTEMR may create NaNs and infinities and
 hence may abort due to a floating point exception in environments
 which do not handle NaNs and infinities in the ieee standard default
 manner.
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.
          For RANGE = 'V' or 'I' and IU - IL < N - 1, SSTEBZ and
          CSTEIN are called
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]A
          A is COMPLEX array, dimension (LDA, N)
          On entry, the Hermitian matrix A.  If UPLO = 'U', the
          leading N-by-N upper triangular part of A contains the
          upper triangular part of the matrix A.  If UPLO = 'L',
          the leading N-by-N lower triangular part of A contains
          the lower triangular part of the matrix A.
          On exit, the lower triangle (if UPLO='L') or the upper
          triangle (if UPLO='U') of A, including the diagonal, is
          destroyed.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[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; 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 REAL
          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 A to tridiagonal form.

          See "Computing Small Singular Values of Bidiagonal Matrices
          with Guaranteed High Relative Accuracy," by Demmel and
          Kahan, LAPACK Working Note #3.

          If high relative accuracy is important, set ABSTOL to
          SLAMCH( 'Safe minimum' ).  Doing so will guarantee that
          eigenvalues are computed to high relative accuracy when
          possible in future releases.  The current code does not
          make any guarantees about high relative accuracy, but
          future releases will. See J. Barlow and J. Demmel,
          "Computing Accurate Eigensystems of Scaled Diagonally
          Dominant Matrices", LAPACK Working Note #7, for a discussion
          of which matrices define their eigenvalues to high relative
          accuracy.
[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', 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).
          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 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]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 eigenvector
          is nonzero only in elements ISUPPZ( 2*i-1 ) through
          ISUPPZ( 2*i ). This is an output of CSTEMR (tridiagonal
          matrix). The support of the eigenvectors of A is typically
          1:N because of the unitary transformations applied by CUNMTR.
          Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1
[out]WORK
          WORK is COMPLEX array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The length of the array WORK.  LWORK >= max(1,2*N).
          For optimal efficiency, LWORK >= (NB+1)*N,
          where NB is the max of the blocksize for CHETRD and for
          CUNMTR as returned by ILAENV.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal 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 REAL array, dimension (MAX(1,LRWORK))
          On exit, if INFO = 0, RWORK(1) returns the optimal
          (and minimal) LRWORK.
[in]LRWORK
          LRWORK is INTEGER
          The length of the array RWORK.  LRWORK >= max(1,24*N).

          If LRWORK = -1, then a workspace query is assumed; the
          routine only calculates the optimal 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 optimal
          (and minimal) LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of the array IWORK.  LIWORK >= max(1,10*N).

          If LIWORK = -1, then a workspace query is assumed; the
          routine only calculates the optimal 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:  Internal error
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Inderjit Dhillon, IBM Almaden, USA
Osni Marques, LBNL/NERSC, USA
Ken Stanley, Computer Science Division, University of California at Berkeley, USA
Jason Riedy, Computer Science Division, University of California at Berkeley, USA

Definition at line 354 of file cheevr.f.

357 *
358 * -- LAPACK driver routine --
359 * -- LAPACK is a software package provided by Univ. of Tennessee, --
360 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
361 *
362 * .. Scalar Arguments ..
363  CHARACTER JOBZ, RANGE, UPLO
364  INTEGER IL, INFO, IU, LDA, LDZ, LIWORK, LRWORK, LWORK,
365  $ M, N
366  REAL ABSTOL, VL, VU
367 * ..
368 * .. Array Arguments ..
369  INTEGER ISUPPZ( * ), IWORK( * )
370  REAL RWORK( * ), W( * )
371  COMPLEX A( LDA, * ), WORK( * ), Z( LDZ, * )
372 * ..
373 *
374 * =====================================================================
375 *
376 * .. Parameters ..
377  REAL ZERO, ONE, TWO
378  parameter( zero = 0.0e+0, one = 1.0e+0, two = 2.0e+0 )
379 * ..
380 * .. Local Scalars ..
381  LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, TEST, VALEIG,
382  $ WANTZ, TRYRAC
383  CHARACTER ORDER
384  INTEGER I, IEEEOK, IINFO, IMAX, INDIBL, INDIFL, INDISP,
385  $ INDIWO, INDRD, INDRDD, INDRE, INDREE, INDRWK,
386  $ INDTAU, INDWK, INDWKN, ISCALE, ITMP1, J, JJ,
387  $ LIWMIN, LLWORK, LLRWORK, LLWRKN, LRWMIN,
388  $ LWKOPT, LWMIN, NB, NSPLIT
389  REAL ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
390  $ SIGMA, SMLNUM, TMP1, VLL, VUU
391 * ..
392 * .. External Functions ..
393  LOGICAL LSAME
394  INTEGER ILAENV
395  REAL CLANSY, SLAMCH
396  EXTERNAL lsame, ilaenv, clansy, slamch
397 * ..
398 * .. External Subroutines ..
399  EXTERNAL chetrd, csscal, cstemr, cstein, cswap, cunmtr,
401 * ..
402 * .. Intrinsic Functions ..
403  INTRINSIC max, min, real, sqrt
404 * ..
405 * .. Executable Statements ..
406 *
407 * Test the input parameters.
408 *
409  ieeeok = ilaenv( 10, 'CHEEVR', 'N', 1, 2, 3, 4 )
410 *
411  lower = lsame( uplo, 'L' )
412  wantz = lsame( jobz, 'V' )
413  alleig = lsame( range, 'A' )
414  valeig = lsame( range, 'V' )
415  indeig = lsame( range, 'I' )
416 *
417  lquery = ( ( lwork.EQ.-1 ) .OR. ( lrwork.EQ.-1 ) .OR.
418  $ ( liwork.EQ.-1 ) )
419 *
420  lrwmin = max( 1, 24*n )
421  liwmin = max( 1, 10*n )
422  lwmin = max( 1, 2*n )
423 *
424  info = 0
425  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
426  info = -1
427  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
428  info = -2
429  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
430  info = -3
431  ELSE IF( n.LT.0 ) THEN
432  info = -4
433  ELSE IF( lda.LT.max( 1, n ) ) THEN
434  info = -6
435  ELSE
436  IF( valeig ) THEN
437  IF( n.GT.0 .AND. vu.LE.vl )
438  $ info = -8
439  ELSE IF( indeig ) THEN
440  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
441  info = -9
442  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
443  info = -10
444  END IF
445  END IF
446  END IF
447  IF( info.EQ.0 ) THEN
448  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
449  info = -15
450  END IF
451  END IF
452 *
453  IF( info.EQ.0 ) THEN
454  nb = ilaenv( 1, 'CHETRD', uplo, n, -1, -1, -1 )
455  nb = max( nb, ilaenv( 1, 'CUNMTR', uplo, n, -1, -1, -1 ) )
456  lwkopt = max( ( nb+1 )*n, lwmin )
457  work( 1 ) = lwkopt
458  rwork( 1 ) = lrwmin
459  iwork( 1 ) = liwmin
460 *
461  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
462  info = -18
463  ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
464  info = -20
465  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
466  info = -22
467  END IF
468  END IF
469 *
470  IF( info.NE.0 ) THEN
471  CALL xerbla( 'CHEEVR', -info )
472  RETURN
473  ELSE IF( lquery ) THEN
474  RETURN
475  END IF
476 *
477 * Quick return if possible
478 *
479  m = 0
480  IF( n.EQ.0 ) THEN
481  work( 1 ) = 1
482  RETURN
483  END IF
484 *
485  IF( n.EQ.1 ) THEN
486  work( 1 ) = 2
487  IF( alleig .OR. indeig ) THEN
488  m = 1
489  w( 1 ) = real( a( 1, 1 ) )
490  ELSE
491  IF( vl.LT.real( a( 1, 1 ) ) .AND. vu.GE.real( a( 1, 1 ) ) )
492  $ THEN
493  m = 1
494  w( 1 ) = real( a( 1, 1 ) )
495  END IF
496  END IF
497  IF( wantz ) THEN
498  z( 1, 1 ) = one
499  isuppz( 1 ) = 1
500  isuppz( 2 ) = 1
501  END IF
502  RETURN
503  END IF
504 *
505 * Get machine constants.
506 *
507  safmin = slamch( 'Safe minimum' )
508  eps = slamch( 'Precision' )
509  smlnum = safmin / eps
510  bignum = one / smlnum
511  rmin = sqrt( smlnum )
512  rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
513 *
514 * Scale matrix to allowable range, if necessary.
515 *
516  iscale = 0
517  abstll = abstol
518  IF (valeig) THEN
519  vll = vl
520  vuu = vu
521  END IF
522  anrm = clansy( 'M', uplo, n, a, lda, rwork )
523  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
524  iscale = 1
525  sigma = rmin / anrm
526  ELSE IF( anrm.GT.rmax ) THEN
527  iscale = 1
528  sigma = rmax / anrm
529  END IF
530  IF( iscale.EQ.1 ) THEN
531  IF( lower ) THEN
532  DO 10 j = 1, n
533  CALL csscal( n-j+1, sigma, a( j, j ), 1 )
534  10 CONTINUE
535  ELSE
536  DO 20 j = 1, n
537  CALL csscal( j, sigma, a( 1, j ), 1 )
538  20 CONTINUE
539  END IF
540  IF( abstol.GT.0 )
541  $ abstll = abstol*sigma
542  IF( valeig ) THEN
543  vll = vl*sigma
544  vuu = vu*sigma
545  END IF
546  END IF
547 
548 * Initialize indices into workspaces. Note: The IWORK indices are
549 * used only if SSTERF or CSTEMR fail.
550 
551 * WORK(INDTAU:INDTAU+N-1) stores the complex scalar factors of the
552 * elementary reflectors used in CHETRD.
553  indtau = 1
554 * INDWK is the starting offset of the remaining complex workspace,
555 * and LLWORK is the remaining complex workspace size.
556  indwk = indtau + n
557  llwork = lwork - indwk + 1
558 
559 * RWORK(INDRD:INDRD+N-1) stores the real tridiagonal's diagonal
560 * entries.
561  indrd = 1
562 * RWORK(INDRE:INDRE+N-1) stores the off-diagonal entries of the
563 * tridiagonal matrix from CHETRD.
564  indre = indrd + n
565 * RWORK(INDRDD:INDRDD+N-1) is a copy of the diagonal entries over
566 * -written by CSTEMR (the SSTERF path copies the diagonal to W).
567  indrdd = indre + n
568 * RWORK(INDREE:INDREE+N-1) is a copy of the off-diagonal entries over
569 * -written while computing the eigenvalues in SSTERF and CSTEMR.
570  indree = indrdd + n
571 * INDRWK is the starting offset of the left-over real workspace, and
572 * LLRWORK is the remaining workspace size.
573  indrwk = indree + n
574  llrwork = lrwork - indrwk + 1
575 
576 * IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in SSTEBZ and
577 * stores the block indices of each of the M<=N eigenvalues.
578  indibl = 1
579 * IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in SSTEBZ and
580 * stores the starting and finishing indices of each block.
581  indisp = indibl + n
582 * IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
583 * that corresponding to eigenvectors that fail to converge in
584 * SSTEIN. This information is discarded; if any fail, the driver
585 * returns INFO > 0.
586  indifl = indisp + n
587 * INDIWO is the offset of the remaining integer workspace.
588  indiwo = indifl + n
589 
590 *
591 * Call CHETRD to reduce Hermitian matrix to tridiagonal form.
592 *
593  CALL chetrd( uplo, n, a, lda, rwork( indrd ), rwork( indre ),
594  $ work( indtau ), work( indwk ), llwork, iinfo )
595 *
596 * If all eigenvalues are desired
597 * then call SSTERF or CSTEMR and CUNMTR.
598 *
599  test = .false.
600  IF( indeig ) THEN
601  IF( il.EQ.1 .AND. iu.EQ.n ) THEN
602  test = .true.
603  END IF
604  END IF
605  IF( ( alleig.OR.test ) .AND. ( ieeeok.EQ.1 ) ) THEN
606  IF( .NOT.wantz ) THEN
607  CALL scopy( n, rwork( indrd ), 1, w, 1 )
608  CALL scopy( n-1, rwork( indre ), 1, rwork( indree ), 1 )
609  CALL ssterf( n, w, rwork( indree ), info )
610  ELSE
611  CALL scopy( n-1, rwork( indre ), 1, rwork( indree ), 1 )
612  CALL scopy( n, rwork( indrd ), 1, rwork( indrdd ), 1 )
613 *
614  IF (abstol .LE. two*n*eps) THEN
615  tryrac = .true.
616  ELSE
617  tryrac = .false.
618  END IF
619  CALL cstemr( jobz, 'A', n, rwork( indrdd ),
620  $ rwork( indree ), vl, vu, il, iu, m, w,
621  $ z, ldz, n, isuppz, tryrac,
622  $ rwork( indrwk ), llrwork,
623  $ iwork, liwork, info )
624 *
625 * Apply unitary matrix used in reduction to tridiagonal
626 * form to eigenvectors returned by CSTEMR.
627 *
628  IF( wantz .AND. info.EQ.0 ) THEN
629  indwkn = indwk
630  llwrkn = lwork - indwkn + 1
631  CALL cunmtr( 'L', uplo, 'N', n, m, a, lda,
632  $ work( indtau ), z, ldz, work( indwkn ),
633  $ llwrkn, iinfo )
634  END IF
635  END IF
636 *
637 *
638  IF( info.EQ.0 ) THEN
639  m = n
640  GO TO 30
641  END IF
642  info = 0
643  END IF
644 *
645 * Otherwise, call SSTEBZ and, if eigenvectors are desired, CSTEIN.
646 * Also call SSTEBZ and CSTEIN if CSTEMR fails.
647 *
648  IF( wantz ) THEN
649  order = 'B'
650  ELSE
651  order = 'E'
652  END IF
653 
654  CALL sstebz( range, order, n, vll, vuu, il, iu, abstll,
655  $ rwork( indrd ), rwork( indre ), m, nsplit, w,
656  $ iwork( indibl ), iwork( indisp ), rwork( indrwk ),
657  $ iwork( indiwo ), info )
658 *
659  IF( wantz ) THEN
660  CALL cstein( n, rwork( indrd ), rwork( indre ), m, w,
661  $ iwork( indibl ), iwork( indisp ), z, ldz,
662  $ rwork( indrwk ), iwork( indiwo ), iwork( indifl ),
663  $ info )
664 *
665 * Apply unitary matrix used in reduction to tridiagonal
666 * form to eigenvectors returned by CSTEIN.
667 *
668  indwkn = indwk
669  llwrkn = lwork - indwkn + 1
670  CALL cunmtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
671  $ ldz, work( indwkn ), llwrkn, iinfo )
672  END IF
673 *
674 * If matrix was scaled, then rescale eigenvalues appropriately.
675 *
676  30 CONTINUE
677  IF( iscale.EQ.1 ) THEN
678  IF( info.EQ.0 ) THEN
679  imax = m
680  ELSE
681  imax = info - 1
682  END IF
683  CALL sscal( imax, one / sigma, w, 1 )
684  END IF
685 *
686 * If eigenvalues are not in order, then sort them, along with
687 * eigenvectors.
688 *
689  IF( wantz ) THEN
690  DO 50 j = 1, m - 1
691  i = 0
692  tmp1 = w( j )
693  DO 40 jj = j + 1, m
694  IF( w( jj ).LT.tmp1 ) THEN
695  i = jj
696  tmp1 = w( jj )
697  END IF
698  40 CONTINUE
699 *
700  IF( i.NE.0 ) THEN
701  itmp1 = iwork( indibl+i-1 )
702  w( i ) = w( j )
703  iwork( indibl+i-1 ) = iwork( indibl+j-1 )
704  w( j ) = tmp1
705  iwork( indibl+j-1 ) = itmp1
706  CALL cswap( n, z( 1, i ), 1, z( 1, j ), 1 )
707  END IF
708  50 CONTINUE
709  END IF
710 *
711 * Set WORK(1) to optimal workspace size.
712 *
713  work( 1 ) = lwkopt
714  rwork( 1 ) = lrwmin
715  iwork( 1 ) = liwmin
716 *
717  RETURN
718 *
719 * End of CHEEVR
720 *
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
ILAENV
Definition: ilaenv.f:162
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:86
subroutine sstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
SSTEBZ
Definition: sstebz.f:273
subroutine csscal(N, SA, CX, INCX)
CSSCAL
Definition: csscal.f:78
subroutine cswap(N, CX, INCX, CY, INCY)
CSWAP
Definition: cswap.f:81
subroutine chetrd(UPLO, N, A, LDA, D, E, TAU, WORK, LWORK, INFO)
CHETRD
Definition: chetrd.f:192
subroutine cunmtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
CUNMTR
Definition: cunmtr.f:172
subroutine cstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
CSTEIN
Definition: cstein.f:182
subroutine cstemr(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, M, W, Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, IWORK, LIWORK, INFO)
CSTEMR
Definition: cstemr.f:338
real function clansy(NORM, UPLO, N, A, LDA, WORK)
CLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition: clansy.f:123
subroutine scopy(N, SX, INCX, SY, INCY)
SCOPY
Definition: scopy.f:82
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:79
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:68
Here is the call graph for this function:
Here is the caller graph for this function: