LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
zheevr.f
Go to the documentation of this file.
1 *> \brief <b> ZHEEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices</b>
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download ZHEEVR + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zheevr.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zheevr.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zheevr.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZHEEVR( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU,
22 * ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK,
23 * RWORK, LRWORK, IWORK, LIWORK, INFO )
24 *
25 * .. Scalar Arguments ..
26 * CHARACTER JOBZ, RANGE, UPLO
27 * INTEGER IL, INFO, IU, LDA, LDZ, LIWORK, LRWORK, LWORK,
28 * $ M, N
29 * DOUBLE PRECISION ABSTOL, VL, VU
30 * ..
31 * .. Array Arguments ..
32 * INTEGER ISUPPZ( * ), IWORK( * )
33 * DOUBLE PRECISION RWORK( * ), W( * )
34 * COMPLEX*16 A( LDA, * ), WORK( * ), Z( LDZ, * )
35 * ..
36 *
37 *
38 *> \par Purpose:
39 * =============
40 *>
41 *> \verbatim
42 *>
43 *> ZHEEVR computes selected eigenvalues and, optionally, eigenvectors
44 *> of a complex Hermitian matrix A. Eigenvalues and eigenvectors can
45 *> be selected by specifying either a range of values or a range of
46 *> indices for the desired eigenvalues.
47 *>
48 *> ZHEEVR first reduces the matrix A to tridiagonal form T with a call
49 *> to ZHETRD. Then, whenever possible, ZHEEVR calls ZSTEMR to compute
50 *> eigenspectrum using Relatively Robust Representations. ZSTEMR
51 *> computes eigenvalues by the dqds algorithm, while orthogonal
52 *> eigenvectors are computed from various "good" L D L^T representations
53 *> (also known as Relatively Robust Representations). Gram-Schmidt
54 *> orthogonalization is avoided as far as possible. More specifically,
55 *> the various steps of the algorithm are as follows.
56 *>
57 *> For each unreduced block (submatrix) of T,
58 *> (a) Compute T - sigma I = L D L^T, so that L and D
59 *> define all the wanted eigenvalues to high relative accuracy.
60 *> This means that small relative changes in the entries of D and L
61 *> cause only small relative changes in the eigenvalues and
62 *> eigenvectors. The standard (unfactored) representation of the
63 *> tridiagonal matrix T does not have this property in general.
64 *> (b) Compute the eigenvalues to suitable accuracy.
65 *> If the eigenvectors are desired, the algorithm attains full
66 *> accuracy of the computed eigenvalues only right before
67 *> the corresponding vectors have to be computed, see steps c) and d).
68 *> (c) For each cluster of close eigenvalues, select a new
69 *> shift close to the cluster, find a new factorization, and refine
70 *> the shifted eigenvalues to suitable accuracy.
71 *> (d) For each eigenvalue with a large enough relative separation compute
72 *> the corresponding eigenvector by forming a rank revealing twisted
73 *> factorization. Go back to (c) for any clusters that remain.
74 *>
75 *> The desired accuracy of the output can be specified by the input
76 *> parameter ABSTOL.
77 *>
78 *> For more details, see DSTEMR's documentation and:
79 *> - Inderjit S. Dhillon and Beresford N. Parlett: "Multiple representations
80 *> to compute orthogonal eigenvectors of symmetric tridiagonal matrices,"
81 *> Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004.
82 *> - Inderjit Dhillon and Beresford Parlett: "Orthogonal Eigenvectors and
83 *> Relative Gaps," SIAM Journal on Matrix Analysis and Applications, Vol. 25,
84 *> 2004. Also LAPACK Working Note 154.
85 *> - Inderjit Dhillon: "A new O(n^2) algorithm for the symmetric
86 *> tridiagonal eigenvalue/eigenvector problem",
87 *> Computer Science Division Technical Report No. UCB/CSD-97-971,
88 *> UC Berkeley, May 1997.
89 *>
90 *>
91 *> Note 1 : ZHEEVR calls ZSTEMR when the full spectrum is requested
92 *> on machines which conform to the ieee-754 floating point standard.
93 *> ZHEEVR calls DSTEBZ and ZSTEIN on non-ieee machines and
94 *> when partial spectrum requests are made.
95 *>
96 *> Normal execution of ZSTEMR may create NaNs and infinities and
97 *> hence may abort due to a floating point exception in environments
98 *> which do not handle NaNs and infinities in the ieee standard default
99 *> manner.
100 *> \endverbatim
101 *
102 * Arguments:
103 * ==========
104 *
105 *> \param[in] JOBZ
106 *> \verbatim
107 *> JOBZ is CHARACTER*1
108 *> = 'N': Compute eigenvalues only;
109 *> = 'V': Compute eigenvalues and eigenvectors.
110 *> \endverbatim
111 *>
112 *> \param[in] RANGE
113 *> \verbatim
114 *> RANGE is CHARACTER*1
115 *> = 'A': all eigenvalues will be found.
116 *> = 'V': all eigenvalues in the half-open interval (VL,VU]
117 *> will be found.
118 *> = 'I': the IL-th through IU-th eigenvalues will be found.
119 *> For RANGE = 'V' or 'I' and IU - IL < N - 1, DSTEBZ and
120 *> ZSTEIN are called
121 *> \endverbatim
122 *>
123 *> \param[in] UPLO
124 *> \verbatim
125 *> UPLO is CHARACTER*1
126 *> = 'U': Upper triangle of A is stored;
127 *> = 'L': Lower triangle of A is stored.
128 *> \endverbatim
129 *>
130 *> \param[in] N
131 *> \verbatim
132 *> N is INTEGER
133 *> The order of the matrix A. N >= 0.
134 *> \endverbatim
135 *>
136 *> \param[in,out] A
137 *> \verbatim
138 *> A is COMPLEX*16 array, dimension (LDA, N)
139 *> On entry, the Hermitian matrix A. If UPLO = 'U', the
140 *> leading N-by-N upper triangular part of A contains the
141 *> upper triangular part of the matrix A. If UPLO = 'L',
142 *> the leading N-by-N lower triangular part of A contains
143 *> the lower triangular part of the matrix A.
144 *> On exit, the lower triangle (if UPLO='L') or the upper
145 *> triangle (if UPLO='U') of A, including the diagonal, is
146 *> destroyed.
147 *> \endverbatim
148 *>
149 *> \param[in] LDA
150 *> \verbatim
151 *> LDA is INTEGER
152 *> The leading dimension of the array A. LDA >= max(1,N).
153 *> \endverbatim
154 *>
155 *> \param[in] VL
156 *> \verbatim
157 *> VL is DOUBLE PRECISION
158 *> If RANGE='V', the lower bound of the interval to
159 *> be searched for eigenvalues. VL < VU.
160 *> Not referenced if RANGE = 'A' or 'I'.
161 *> \endverbatim
162 *>
163 *> \param[in] VU
164 *> \verbatim
165 *> VU is DOUBLE PRECISION
166 *> If RANGE='V', the upper bound of the interval to
167 *> be searched for eigenvalues. VL < VU.
168 *> Not referenced if RANGE = 'A' or 'I'.
169 *> \endverbatim
170 *>
171 *> \param[in] IL
172 *> \verbatim
173 *> IL is INTEGER
174 *> If RANGE='I', the index of the
175 *> smallest eigenvalue to be returned.
176 *> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
177 *> Not referenced if RANGE = 'A' or 'V'.
178 *> \endverbatim
179 *>
180 *> \param[in] IU
181 *> \verbatim
182 *> IU is INTEGER
183 *> If RANGE='I', the index of the
184 *> largest eigenvalue to be returned.
185 *> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
186 *> Not referenced if RANGE = 'A' or 'V'.
187 *> \endverbatim
188 *>
189 *> \param[in] ABSTOL
190 *> \verbatim
191 *> ABSTOL is DOUBLE PRECISION
192 *> The absolute error tolerance for the eigenvalues.
193 *> An approximate eigenvalue is accepted as converged
194 *> when it is determined to lie in an interval [a,b]
195 *> of width less than or equal to
196 *>
197 *> ABSTOL + EPS * max( |a|,|b| ) ,
198 *>
199 *> where EPS is the machine precision. If ABSTOL is less than
200 *> or equal to zero, then EPS*|T| will be used in its place,
201 *> where |T| is the 1-norm of the tridiagonal matrix obtained
202 *> by reducing A to tridiagonal form.
203 *>
204 *> See "Computing Small Singular Values of Bidiagonal Matrices
205 *> with Guaranteed High Relative Accuracy," by Demmel and
206 *> Kahan, LAPACK Working Note #3.
207 *>
208 *> If high relative accuracy is important, set ABSTOL to
209 *> DLAMCH( 'Safe minimum' ). Doing so will guarantee that
210 *> eigenvalues are computed to high relative accuracy when
211 *> possible in future releases. The current code does not
212 *> make any guarantees about high relative accuracy, but
213 *> furutre releases will. See J. Barlow and J. Demmel,
214 *> "Computing Accurate Eigensystems of Scaled Diagonally
215 *> Dominant Matrices", LAPACK Working Note #7, for a discussion
216 *> of which matrices define their eigenvalues to high relative
217 *> accuracy.
218 *> \endverbatim
219 *>
220 *> \param[out] M
221 *> \verbatim
222 *> M is INTEGER
223 *> The total number of eigenvalues found. 0 <= M <= N.
224 *> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
225 *> \endverbatim
226 *>
227 *> \param[out] W
228 *> \verbatim
229 *> W is DOUBLE PRECISION array, dimension (N)
230 *> The first M elements contain the selected eigenvalues in
231 *> ascending order.
232 *> \endverbatim
233 *>
234 *> \param[out] Z
235 *> \verbatim
236 *> Z is COMPLEX*16 array, dimension (LDZ, max(1,M))
237 *> If JOBZ = 'V', then if INFO = 0, the first M columns of Z
238 *> contain the orthonormal eigenvectors of the matrix A
239 *> corresponding to the selected eigenvalues, with the i-th
240 *> column of Z holding the eigenvector associated with W(i).
241 *> If JOBZ = 'N', then Z is not referenced.
242 *> Note: the user must ensure that at least max(1,M) columns are
243 *> supplied in the array Z; if RANGE = 'V', the exact value of M
244 *> is not known in advance and an upper bound must be used.
245 *> \endverbatim
246 *>
247 *> \param[in] LDZ
248 *> \verbatim
249 *> LDZ is INTEGER
250 *> The leading dimension of the array Z. LDZ >= 1, and if
251 *> JOBZ = 'V', LDZ >= max(1,N).
252 *> \endverbatim
253 *>
254 *> \param[out] ISUPPZ
255 *> \verbatim
256 *> ISUPPZ is INTEGER array, dimension ( 2*max(1,M) )
257 *> The support of the eigenvectors in Z, i.e., the indices
258 *> indicating the nonzero elements in Z. The i-th eigenvector
259 *> is nonzero only in elements ISUPPZ( 2*i-1 ) through
260 *> ISUPPZ( 2*i ).
261 *> Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1
262 *> \endverbatim
263 *>
264 *> \param[out] WORK
265 *> \verbatim
266 *> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
267 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
268 *> \endverbatim
269 *>
270 *> \param[in] LWORK
271 *> \verbatim
272 *> LWORK is INTEGER
273 *> The length of the array WORK. LWORK >= max(1,2*N).
274 *> For optimal efficiency, LWORK >= (NB+1)*N,
275 *> where NB is the max of the blocksize for ZHETRD and for
276 *> ZUNMTR as returned by ILAENV.
277 *>
278 *> If LWORK = -1, then a workspace query is assumed; the routine
279 *> only calculates the optimal sizes of the WORK, RWORK and
280 *> IWORK arrays, returns these values as the first entries of
281 *> the WORK, RWORK and IWORK arrays, and no error message
282 *> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
283 *> \endverbatim
284 *>
285 *> \param[out] RWORK
286 *> \verbatim
287 *> RWORK is DOUBLE PRECISION array, dimension (MAX(1,LRWORK))
288 *> On exit, if INFO = 0, RWORK(1) returns the optimal
289 *> (and minimal) LRWORK.
290 *> \endverbatim
291 *>
292 *> \param[in] LRWORK
293 *> \verbatim
294 *> LRWORK is INTEGER
295 *> The length of the array RWORK. LRWORK >= max(1,24*N).
296 *>
297 *> If LRWORK = -1, then a workspace query is assumed; the
298 *> routine only calculates the optimal sizes of the WORK, RWORK
299 *> and IWORK arrays, returns these values as the first entries
300 *> of the WORK, RWORK and IWORK arrays, and no error message
301 *> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
302 *> \endverbatim
303 *>
304 *> \param[out] IWORK
305 *> \verbatim
306 *> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
307 *> On exit, if INFO = 0, IWORK(1) returns the optimal
308 *> (and minimal) LIWORK.
309 *> \endverbatim
310 *>
311 *> \param[in] LIWORK
312 *> \verbatim
313 *> LIWORK is INTEGER
314 *> The dimension of the array IWORK. LIWORK >= max(1,10*N).
315 *>
316 *> If LIWORK = -1, then a workspace query is assumed; the
317 *> routine only calculates the optimal sizes of the WORK, RWORK
318 *> and IWORK arrays, returns these values as the first entries
319 *> of the WORK, RWORK and IWORK arrays, and no error message
320 *> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
321 *> \endverbatim
322 *>
323 *> \param[out] INFO
324 *> \verbatim
325 *> INFO is INTEGER
326 *> = 0: successful exit
327 *> < 0: if INFO = -i, the i-th argument had an illegal value
328 *> > 0: Internal error
329 *> \endverbatim
330 *
331 * Authors:
332 * ========
333 *
334 *> \author Univ. of Tennessee
335 *> \author Univ. of California Berkeley
336 *> \author Univ. of Colorado Denver
337 *> \author NAG Ltd.
338 *
339 *> \date June 2016
340 *
341 *> \ingroup complex16HEeigen
342 *
343 *> \par Contributors:
344 * ==================
345 *>
346 *> Inderjit Dhillon, IBM Almaden, USA \n
347 *> Osni Marques, LBNL/NERSC, USA \n
348 *> Ken Stanley, Computer Science Division, University of
349 *> California at Berkeley, USA \n
350 *> Jason Riedy, Computer Science Division, University of
351 *> California at Berkeley, USA \n
352 *>
353 * =====================================================================
354  SUBROUTINE zheevr( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU,
355  $ abstol, m, w, z, ldz, isuppz, work, lwork,
356  $ rwork, lrwork, iwork, liwork, info )
357 *
358 * -- LAPACK driver routine (version 3.6.1) --
359 * -- LAPACK is a software package provided by Univ. of Tennessee, --
360 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
361 * June 2016
362 *
363 * .. Scalar Arguments ..
364  CHARACTER JOBZ, RANGE, UPLO
365  INTEGER IL, INFO, IU, LDA, LDZ, LIWORK, LRWORK, LWORK,
366  $ m, n
367  DOUBLE PRECISION ABSTOL, VL, VU
368 * ..
369 * .. Array Arguments ..
370  INTEGER ISUPPZ( * ), IWORK( * )
371  DOUBLE PRECISION RWORK( * ), W( * )
372  COMPLEX*16 A( lda, * ), WORK( * ), Z( ldz, * )
373 * ..
374 *
375 * =====================================================================
376 *
377 * .. Parameters ..
378  DOUBLE PRECISION ZERO, ONE, TWO
379  parameter ( zero = 0.0d+0, one = 1.0d+0, two = 2.0d+0 )
380 * ..
381 * .. Local Scalars ..
382  LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, TEST, VALEIG,
383  $ wantz, tryrac
384  CHARACTER ORDER
385  INTEGER I, IEEEOK, IINFO, IMAX, INDIBL, INDIFL, INDISP,
386  $ indiwo, indrd, indrdd, indre, indree, indrwk,
387  $ indtau, indwk, indwkn, iscale, itmp1, j, jj,
388  $ liwmin, llwork, llrwork, llwrkn, lrwmin,
389  $ lwkopt, lwmin, nb, nsplit
390  DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
391  $ sigma, smlnum, tmp1, vll, vuu
392 * ..
393 * .. External Functions ..
394  LOGICAL LSAME
395  INTEGER ILAENV
396  DOUBLE PRECISION DLAMCH, ZLANSY
397  EXTERNAL lsame, ilaenv, dlamch, zlansy
398 * ..
399 * .. External Subroutines ..
400  EXTERNAL dcopy, dscal, dstebz, dsterf, xerbla, zdscal,
402 * ..
403 * .. Intrinsic Functions ..
404  INTRINSIC dble, max, min, sqrt
405 * ..
406 * .. Executable Statements ..
407 *
408 * Test the input parameters.
409 *
410  ieeeok = ilaenv( 10, 'ZHEEVR', 'N', 1, 2, 3, 4 )
411 *
412  lower = lsame( uplo, 'L' )
413  wantz = lsame( jobz, 'V' )
414  alleig = lsame( range, 'A' )
415  valeig = lsame( range, 'V' )
416  indeig = lsame( range, 'I' )
417 *
418  lquery = ( ( lwork.EQ.-1 ) .OR. ( lrwork.EQ.-1 ) .OR.
419  $ ( liwork.EQ.-1 ) )
420 *
421  lrwmin = max( 1, 24*n )
422  liwmin = max( 1, 10*n )
423  lwmin = max( 1, 2*n )
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( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
431  info = -3
432  ELSE IF( n.LT.0 ) THEN
433  info = -4
434  ELSE IF( lda.LT.max( 1, n ) ) THEN
435  info = -6
436  ELSE
437  IF( valeig ) THEN
438  IF( n.GT.0 .AND. vu.LE.vl )
439  $ info = -8
440  ELSE IF( indeig ) THEN
441  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
442  info = -9
443  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
444  info = -10
445  END IF
446  END IF
447  END IF
448  IF( info.EQ.0 ) THEN
449  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
450  info = -15
451  END IF
452  END IF
453 *
454  IF( info.EQ.0 ) THEN
455  nb = ilaenv( 1, 'ZHETRD', uplo, n, -1, -1, -1 )
456  nb = max( nb, ilaenv( 1, 'ZUNMTR', uplo, n, -1, -1, -1 ) )
457  lwkopt = max( ( nb+1 )*n, lwmin )
458  work( 1 ) = lwkopt
459  rwork( 1 ) = lrwmin
460  iwork( 1 ) = liwmin
461 *
462  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
463  info = -18
464  ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
465  info = -20
466  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
467  info = -22
468  END IF
469  END IF
470 *
471  IF( info.NE.0 ) THEN
472  CALL xerbla( 'ZHEEVR', -info )
473  RETURN
474  ELSE IF( lquery ) THEN
475  RETURN
476  END IF
477 *
478 * Quick return if possible
479 *
480  m = 0
481  IF( n.EQ.0 ) THEN
482  work( 1 ) = 1
483  RETURN
484  END IF
485 *
486  IF( n.EQ.1 ) THEN
487  work( 1 ) = 2
488  IF( alleig .OR. indeig ) THEN
489  m = 1
490  w( 1 ) = dble( a( 1, 1 ) )
491  ELSE
492  IF( vl.LT.dble( a( 1, 1 ) ) .AND. vu.GE.dble( a( 1, 1 ) ) )
493  $ THEN
494  m = 1
495  w( 1 ) = dble( a( 1, 1 ) )
496  END IF
497  END IF
498  IF( wantz ) THEN
499  z( 1, 1 ) = one
500  isuppz( 1 ) = 1
501  isuppz( 2 ) = 1
502  END IF
503  RETURN
504  END IF
505 *
506 * Get machine constants.
507 *
508  safmin = dlamch( 'Safe minimum' )
509  eps = dlamch( 'Precision' )
510  smlnum = safmin / eps
511  bignum = one / smlnum
512  rmin = sqrt( smlnum )
513  rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
514 *
515 * Scale matrix to allowable range, if necessary.
516 *
517  iscale = 0
518  abstll = abstol
519  IF (valeig) THEN
520  vll = vl
521  vuu = vu
522  END IF
523  anrm = zlansy( 'M', uplo, n, a, lda, rwork )
524  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
525  iscale = 1
526  sigma = rmin / anrm
527  ELSE IF( anrm.GT.rmax ) THEN
528  iscale = 1
529  sigma = rmax / anrm
530  END IF
531  IF( iscale.EQ.1 ) THEN
532  IF( lower ) THEN
533  DO 10 j = 1, n
534  CALL zdscal( n-j+1, sigma, a( j, j ), 1 )
535  10 CONTINUE
536  ELSE
537  DO 20 j = 1, n
538  CALL zdscal( j, sigma, a( 1, j ), 1 )
539  20 CONTINUE
540  END IF
541  IF( abstol.GT.0 )
542  $ abstll = abstol*sigma
543  IF( valeig ) THEN
544  vll = vl*sigma
545  vuu = vu*sigma
546  END IF
547  END IF
548 
549 * Initialize indices into workspaces. Note: The IWORK indices are
550 * used only if DSTERF or ZSTEMR fail.
551 
552 * WORK(INDTAU:INDTAU+N-1) stores the complex scalar factors of the
553 * elementary reflectors used in ZHETRD.
554  indtau = 1
555 * INDWK is the starting offset of the remaining complex workspace,
556 * and LLWORK is the remaining complex workspace size.
557  indwk = indtau + n
558  llwork = lwork - indwk + 1
559 
560 * RWORK(INDRD:INDRD+N-1) stores the real tridiagonal's diagonal
561 * entries.
562  indrd = 1
563 * RWORK(INDRE:INDRE+N-1) stores the off-diagonal entries of the
564 * tridiagonal matrix from ZHETRD.
565  indre = indrd + n
566 * RWORK(INDRDD:INDRDD+N-1) is a copy of the diagonal entries over
567 * -written by ZSTEMR (the DSTERF path copies the diagonal to W).
568  indrdd = indre + n
569 * RWORK(INDREE:INDREE+N-1) is a copy of the off-diagonal entries over
570 * -written while computing the eigenvalues in DSTERF and ZSTEMR.
571  indree = indrdd + n
572 * INDRWK is the starting offset of the left-over real workspace, and
573 * LLRWORK is the remaining workspace size.
574  indrwk = indree + n
575  llrwork = lrwork - indrwk + 1
576 
577 * IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in DSTEBZ and
578 * stores the block indices of each of the M<=N eigenvalues.
579  indibl = 1
580 * IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in DSTEBZ and
581 * stores the starting and finishing indices of each block.
582  indisp = indibl + n
583 * IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
584 * that corresponding to eigenvectors that fail to converge in
585 * DSTEIN. This information is discarded; if any fail, the driver
586 * returns INFO > 0.
587  indifl = indisp + n
588 * INDIWO is the offset of the remaining integer workspace.
589  indiwo = indifl + n
590 
591 *
592 * Call ZHETRD to reduce Hermitian matrix to tridiagonal form.
593 *
594  CALL zhetrd( uplo, n, a, lda, rwork( indrd ), rwork( indre ),
595  $ work( indtau ), work( indwk ), llwork, iinfo )
596 *
597 * If all eigenvalues are desired
598 * then call DSTERF or ZSTEMR and ZUNMTR.
599 *
600  test = .false.
601  IF( indeig ) THEN
602  IF( il.EQ.1 .AND. iu.EQ.n ) THEN
603  test = .true.
604  END IF
605  END IF
606  IF( ( alleig.OR.test ) .AND. ( ieeeok.EQ.1 ) ) THEN
607  IF( .NOT.wantz ) THEN
608  CALL dcopy( n, rwork( indrd ), 1, w, 1 )
609  CALL dcopy( n-1, rwork( indre ), 1, rwork( indree ), 1 )
610  CALL dsterf( n, w, rwork( indree ), info )
611  ELSE
612  CALL dcopy( n-1, rwork( indre ), 1, rwork( indree ), 1 )
613  CALL dcopy( n, rwork( indrd ), 1, rwork( indrdd ), 1 )
614 *
615  IF (abstol .LE. two*n*eps) THEN
616  tryrac = .true.
617  ELSE
618  tryrac = .false.
619  END IF
620  CALL zstemr( jobz, 'A', n, rwork( indrdd ),
621  $ rwork( indree ), vl, vu, il, iu, m, w,
622  $ z, ldz, n, isuppz, tryrac,
623  $ rwork( indrwk ), llrwork,
624  $ iwork, liwork, info )
625 *
626 * Apply unitary matrix used in reduction to tridiagonal
627 * form to eigenvectors returned by ZSTEIN.
628 *
629  IF( wantz .AND. info.EQ.0 ) THEN
630  indwkn = indwk
631  llwrkn = lwork - indwkn + 1
632  CALL zunmtr( 'L', uplo, 'N', n, m, a, lda,
633  $ work( indtau ), z, ldz, work( indwkn ),
634  $ llwrkn, iinfo )
635  END IF
636  END IF
637 *
638 *
639  IF( info.EQ.0 ) THEN
640  m = n
641  GO TO 30
642  END IF
643  info = 0
644  END IF
645 *
646 * Otherwise, call DSTEBZ and, if eigenvectors are desired, ZSTEIN.
647 * Also call DSTEBZ and ZSTEIN if ZSTEMR fails.
648 *
649  IF( wantz ) THEN
650  order = 'B'
651  ELSE
652  order = 'E'
653  END IF
654 
655  CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
656  $ rwork( indrd ), rwork( indre ), m, nsplit, w,
657  $ iwork( indibl ), iwork( indisp ), rwork( indrwk ),
658  $ iwork( indiwo ), info )
659 *
660  IF( wantz ) THEN
661  CALL zstein( n, rwork( indrd ), rwork( indre ), m, w,
662  $ iwork( indibl ), iwork( indisp ), z, ldz,
663  $ rwork( indrwk ), iwork( indiwo ), iwork( indifl ),
664  $ info )
665 *
666 * Apply unitary matrix used in reduction to tridiagonal
667 * form to eigenvectors returned by ZSTEIN.
668 *
669  indwkn = indwk
670  llwrkn = lwork - indwkn + 1
671  CALL zunmtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
672  $ ldz, work( indwkn ), llwrkn, iinfo )
673  END IF
674 *
675 * If matrix was scaled, then rescale eigenvalues appropriately.
676 *
677  30 CONTINUE
678  IF( iscale.EQ.1 ) THEN
679  IF( info.EQ.0 ) THEN
680  imax = m
681  ELSE
682  imax = info - 1
683  END IF
684  CALL dscal( imax, one / sigma, w, 1 )
685  END IF
686 *
687 * If eigenvalues are not in order, then sort them, along with
688 * eigenvectors.
689 *
690  IF( wantz ) THEN
691  DO 50 j = 1, m - 1
692  i = 0
693  tmp1 = w( j )
694  DO 40 jj = j + 1, m
695  IF( w( jj ).LT.tmp1 ) THEN
696  i = jj
697  tmp1 = w( jj )
698  END IF
699  40 CONTINUE
700 *
701  IF( i.NE.0 ) THEN
702  itmp1 = iwork( indibl+i-1 )
703  w( i ) = w( j )
704  iwork( indibl+i-1 ) = iwork( indibl+j-1 )
705  w( j ) = tmp1
706  iwork( indibl+j-1 ) = itmp1
707  CALL zswap( n, z( 1, i ), 1, z( 1, j ), 1 )
708  END IF
709  50 CONTINUE
710  END IF
711 *
712 * Set WORK(1) to optimal workspace size.
713 *
714  work( 1 ) = lwkopt
715  rwork( 1 ) = lrwmin
716  iwork( 1 ) = liwmin
717 *
718  RETURN
719 *
720 * End of ZHEEVR
721 *
722  END
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
subroutine dstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
DSTEBZ
Definition: dstebz.f:275
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine zhetrd(UPLO, N, A, LDA, D, E, TAU, WORK, LWORK, INFO)
ZHETRD
Definition: zhetrd.f:194
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:52
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
ZSTEIN
Definition: zstein.f:184
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
subroutine zheevr(JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
ZHEEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices ...
Definition: zheevr.f:357
subroutine zunmtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
ZUNMTR
Definition: zunmtr.f:173
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:54
subroutine zstemr(JOBZ, RANGE, N, D, E, VL, VU, IL, IU, M, W, Z, LDZ, NZC, ISUPPZ, TRYRAC, WORK, LWORK, IWORK, LIWORK, INFO)
ZSTEMR
Definition: zstemr.f:340