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