LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
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 *> \endverbatim
157 *>
158 *> \param[in] VU
159 *> \verbatim
160 *> VU is DOUBLE PRECISION
161 *> If RANGE='V', the lower and upper bounds of the interval to
162 *> be searched for eigenvalues. VL < VU.
163 *> Not referenced if RANGE = 'A' or 'I'.
164 *> \endverbatim
165 *>
166 *> \param[in] IL
167 *> \verbatim
168 *> IL is INTEGER
169 *> \endverbatim
170 *>
171 *> \param[in] IU
172 *> \verbatim
173 *> IU is INTEGER
174 *> If RANGE='I', the indices (in ascending order) of the
175 *> smallest and largest eigenvalues 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] ABSTOL
181 *> \verbatim
182 *> ABSTOL is DOUBLE PRECISION
183 *> The absolute error tolerance for the eigenvalues.
184 *> An approximate eigenvalue is accepted as converged
185 *> when it is determined to lie in an interval [a,b]
186 *> of width less than or equal to
187 *>
188 *> ABSTOL + EPS * max( |a|,|b| ) ,
189 *>
190 *> where EPS is the machine precision. If ABSTOL is less than
191 *> or equal to zero, then EPS*|T| will be used in its place,
192 *> where |T| is the 1-norm of the tridiagonal matrix obtained
193 *> by reducing A to tridiagonal form.
194 *>
195 *> See "Computing Small Singular Values of Bidiagonal Matrices
196 *> with Guaranteed High Relative Accuracy," by Demmel and
197 *> Kahan, LAPACK Working Note #3.
198 *>
199 *> If high relative accuracy is important, set ABSTOL to
200 *> DLAMCH( 'Safe minimum' ). Doing so will guarantee that
201 *> eigenvalues are computed to high relative accuracy when
202 *> possible in future releases. The current code does not
203 *> make any guarantees about high relative accuracy, but
204 *> future releases will. See J. Barlow and J. Demmel,
205 *> "Computing Accurate Eigensystems of Scaled Diagonally
206 *> Dominant Matrices", LAPACK Working Note #7, for a discussion
207 *> of which matrices define their eigenvalues to high relative
208 *> accuracy.
209 *> \endverbatim
210 *>
211 *> \param[out] M
212 *> \verbatim
213 *> M is INTEGER
214 *> The total number of eigenvalues found. 0 <= M <= N.
215 *> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
216 *> \endverbatim
217 *>
218 *> \param[out] W
219 *> \verbatim
220 *> W is DOUBLE PRECISION array, dimension (N)
221 *> The first M elements contain the selected eigenvalues in
222 *> ascending order.
223 *> \endverbatim
224 *>
225 *> \param[out] Z
226 *> \verbatim
227 *> Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M))
228 *> If JOBZ = 'V', then if INFO = 0, the first M columns of Z
229 *> contain the orthonormal eigenvectors of the matrix A
230 *> corresponding to the selected eigenvalues, with the i-th
231 *> column of Z holding the eigenvector associated with W(i).
232 *> If JOBZ = 'N', then Z is not referenced.
233 *> Note: the user must ensure that at least max(1,M) columns are
234 *> supplied in the array Z; if RANGE = 'V', the exact value of M
235 *> is not known in advance and an upper bound must be used.
236 *> Supplying N columns is always safe.
237 *> \endverbatim
238 *>
239 *> \param[in] LDZ
240 *> \verbatim
241 *> LDZ is INTEGER
242 *> The leading dimension of the array Z. LDZ >= 1, and if
243 *> JOBZ = 'V', LDZ >= max(1,N).
244 *> \endverbatim
245 *>
246 *> \param[out] ISUPPZ
247 *> \verbatim
248 *> ISUPPZ is INTEGER array, dimension ( 2*max(1,M) )
249 *> The support of the eigenvectors in Z, i.e., the indices
250 *> indicating the nonzero elements in Z. The i-th eigenvector
251 *> is nonzero only in elements ISUPPZ( 2*i-1 ) through
252 *> ISUPPZ( 2*i ).
253 *> Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1
254 *> \endverbatim
255 *>
256 *> \param[out] WORK
257 *> \verbatim
258 *> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
259 *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
260 *> \endverbatim
261 *>
262 *> \param[in] LWORK
263 *> \verbatim
264 *> LWORK is INTEGER
265 *> The dimension of the array WORK. LWORK >= max(1,26*N).
266 *> For optimal efficiency, LWORK >= (NB+6)*N,
267 *> where NB is the max of the blocksize for DSYTRD and DORMTR
268 *> returned by ILAENV.
269 *>
270 *> If LWORK = -1, then a workspace query is assumed; the routine
271 *> only calculates the optimal size of the WORK array, returns
272 *> this value as the first entry of the WORK array, and no error
273 *> message related to LWORK is issued by XERBLA.
274 *> \endverbatim
275 *>
276 *> \param[out] IWORK
277 *> \verbatim
278 *> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
279 *> On exit, if INFO = 0, IWORK(1) returns the optimal LWORK.
280 *> \endverbatim
281 *>
282 *> \param[in] LIWORK
283 *> \verbatim
284 *> LIWORK is INTEGER
285 *> The dimension of the array IWORK. LIWORK >= max(1,10*N).
286 *>
287 *> If LIWORK = -1, then a workspace query is assumed; the
288 *> routine only calculates the optimal size of the IWORK array,
289 *> returns this value as the first entry of the IWORK array, and
290 *> no error message related to LIWORK is issued by XERBLA.
291 *> \endverbatim
292 *>
293 *> \param[out] INFO
294 *> \verbatim
295 *> INFO is INTEGER
296 *> = 0: successful exit
297 *> < 0: if INFO = -i, the i-th argument had an illegal value
298 *> > 0: Internal error
299 *> \endverbatim
300 *
301 * Authors:
302 * ========
303 *
304 *> \author Univ. of Tennessee
305 *> \author Univ. of California Berkeley
306 *> \author Univ. of Colorado Denver
307 *> \author NAG Ltd.
308 *
309 *> \date September 2012
310 *
311 *> \ingroup doubleSYeigen
312 *
313 *> \par Contributors:
314 * ==================
315 *>
316 *> Inderjit Dhillon, IBM Almaden, USA \n
317 *> Osni Marques, LBNL/NERSC, USA \n
318 *> Ken Stanley, Computer Science Division, University of
319 *> California at Berkeley, USA \n
320 *> Jason Riedy, Computer Science Division, University of
321 *> California at Berkeley, USA \n
322 *>
323 * =====================================================================
324  SUBROUTINE dsyevr( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU,
325  $ abstol, m, w, z, ldz, isuppz, work, lwork,
326  $ iwork, liwork, info )
327 *
328 * -- LAPACK driver routine (version 3.4.2) --
329 * -- LAPACK is a software package provided by Univ. of Tennessee, --
330 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
331 * September 2012
332 *
333 * .. Scalar Arguments ..
334  CHARACTER jobz, range, uplo
335  INTEGER il, info, iu, lda, ldz, liwork, lwork, m, n
336  DOUBLE PRECISION abstol, vl, vu
337 * ..
338 * .. Array Arguments ..
339  INTEGER isuppz( * ), iwork( * )
340  DOUBLE PRECISION a( lda, * ), w( * ), work( * ), z( ldz, * )
341 * ..
342 *
343 * =====================================================================
344 *
345 * .. Parameters ..
346  DOUBLE PRECISION zero, one, two
347  parameter( zero = 0.0d+0, one = 1.0d+0, two = 2.0d+0 )
348 * ..
349 * .. Local Scalars ..
350  LOGICAL alleig, indeig, lower, lquery, valeig, wantz,
351  $ tryrac
352  CHARACTER order
353  INTEGER i, ieeeok, iinfo, imax, indd, inddd, inde,
354  $ indee, indibl, indifl, indisp, indiwo, indtau,
355  $ indwk, indwkn, iscale, j, jj, liwmin,
356  $ llwork, llwrkn, lwkopt, lwmin, nb, nsplit
357  DOUBLE PRECISION abstll, anrm, bignum, eps, rmax, rmin, safmin,
358  $ sigma, smlnum, tmp1, vll, vuu
359 * ..
360 * .. External Functions ..
361  LOGICAL lsame
362  INTEGER ilaenv
363  DOUBLE PRECISION dlamch, dlansy
364  EXTERNAL lsame, ilaenv, dlamch, dlansy
365 * ..
366 * .. External Subroutines ..
367  EXTERNAL dcopy, dormtr, dscal, dstebz, dstemr, dstein,
369 * ..
370 * .. Intrinsic Functions ..
371  INTRINSIC max, min, sqrt
372 * ..
373 * .. Executable Statements ..
374 *
375 * Test the input parameters.
376 *
377  ieeeok = ilaenv( 10, 'DSYEVR', 'N', 1, 2, 3, 4 )
378 *
379  lower = lsame( uplo, 'L' )
380  wantz = lsame( jobz, 'V' )
381  alleig = lsame( range, 'A' )
382  valeig = lsame( range, 'V' )
383  indeig = lsame( range, 'I' )
384 *
385  lquery = ( ( lwork.EQ.-1 ) .OR. ( liwork.EQ.-1 ) )
386 *
387  lwmin = max( 1, 26*n )
388  liwmin = max( 1, 10*n )
389 *
390  info = 0
391  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
392  info = -1
393  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
394  info = -2
395  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
396  info = -3
397  ELSE IF( n.LT.0 ) THEN
398  info = -4
399  ELSE IF( lda.LT.max( 1, n ) ) THEN
400  info = -6
401  ELSE
402  IF( valeig ) THEN
403  IF( n.GT.0 .AND. vu.LE.vl )
404  $ info = -8
405  ELSE IF( indeig ) THEN
406  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
407  info = -9
408  ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
409  info = -10
410  END IF
411  END IF
412  END IF
413  IF( info.EQ.0 ) THEN
414  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
415  info = -15
416  ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
417  info = -18
418  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
419  info = -20
420  END IF
421  END IF
422 *
423  IF( info.EQ.0 ) THEN
424  nb = ilaenv( 1, 'DSYTRD', uplo, n, -1, -1, -1 )
425  nb = max( nb, ilaenv( 1, 'DORMTR', uplo, n, -1, -1, -1 ) )
426  lwkopt = max( ( nb+1 )*n, lwmin )
427  work( 1 ) = lwkopt
428  iwork( 1 ) = liwmin
429  END IF
430 *
431  IF( info.NE.0 ) THEN
432  CALL xerbla( 'DSYEVR', -info )
433  return
434  ELSE IF( lquery ) THEN
435  return
436  END IF
437 *
438 * Quick return if possible
439 *
440  m = 0
441  IF( n.EQ.0 ) THEN
442  work( 1 ) = 1
443  return
444  END IF
445 *
446  IF( n.EQ.1 ) THEN
447  work( 1 ) = 7
448  IF( alleig .OR. indeig ) THEN
449  m = 1
450  w( 1 ) = a( 1, 1 )
451  ELSE
452  IF( vl.LT.a( 1, 1 ) .AND. vu.GE.a( 1, 1 ) ) THEN
453  m = 1
454  w( 1 ) = a( 1, 1 )
455  END IF
456  END IF
457  IF( wantz ) THEN
458  z( 1, 1 ) = one
459  isuppz( 1 ) = 1
460  isuppz( 2 ) = 1
461  END IF
462  return
463  END IF
464 *
465 * Get machine constants.
466 *
467  safmin = dlamch( 'Safe minimum' )
468  eps = dlamch( 'Precision' )
469  smlnum = safmin / eps
470  bignum = one / smlnum
471  rmin = sqrt( smlnum )
472  rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
473 *
474 * Scale matrix to allowable range, if necessary.
475 *
476  iscale = 0
477  abstll = abstol
478  IF (valeig) THEN
479  vll = vl
480  vuu = vu
481  END IF
482  anrm = dlansy( 'M', uplo, n, a, lda, work )
483  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
484  iscale = 1
485  sigma = rmin / anrm
486  ELSE IF( anrm.GT.rmax ) THEN
487  iscale = 1
488  sigma = rmax / anrm
489  END IF
490  IF( iscale.EQ.1 ) THEN
491  IF( lower ) THEN
492  DO 10 j = 1, n
493  CALL dscal( n-j+1, sigma, a( j, j ), 1 )
494  10 continue
495  ELSE
496  DO 20 j = 1, n
497  CALL dscal( j, sigma, a( 1, j ), 1 )
498  20 continue
499  END IF
500  IF( abstol.GT.0 )
501  $ abstll = abstol*sigma
502  IF( valeig ) THEN
503  vll = vl*sigma
504  vuu = vu*sigma
505  END IF
506  END IF
507 
508 * Initialize indices into workspaces. Note: The IWORK indices are
509 * used only if DSTERF or DSTEMR fail.
510 
511 * WORK(INDTAU:INDTAU+N-1) stores the scalar factors of the
512 * elementary reflectors used in DSYTRD.
513  indtau = 1
514 * WORK(INDD:INDD+N-1) stores the tridiagonal's diagonal entries.
515  indd = indtau + n
516 * WORK(INDE:INDE+N-1) stores the off-diagonal entries of the
517 * tridiagonal matrix from DSYTRD.
518  inde = indd + n
519 * WORK(INDDD:INDDD+N-1) is a copy of the diagonal entries over
520 * -written by DSTEMR (the DSTERF path copies the diagonal to W).
521  inddd = inde + n
522 * WORK(INDEE:INDEE+N-1) is a copy of the off-diagonal entries over
523 * -written while computing the eigenvalues in DSTERF and DSTEMR.
524  indee = inddd + n
525 * INDWK is the starting offset of the left-over workspace, and
526 * LLWORK is the remaining workspace size.
527  indwk = indee + n
528  llwork = lwork - indwk + 1
529 
530 * IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in DSTEBZ and
531 * stores the block indices of each of the M<=N eigenvalues.
532  indibl = 1
533 * IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in DSTEBZ and
534 * stores the starting and finishing indices of each block.
535  indisp = indibl + n
536 * IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
537 * that corresponding to eigenvectors that fail to converge in
538 * DSTEIN. This information is discarded; if any fail, the driver
539 * returns INFO > 0.
540  indifl = indisp + n
541 * INDIWO is the offset of the remaining integer workspace.
542  indiwo = indifl + n
543 
544 *
545 * Call DSYTRD to reduce symmetric matrix to tridiagonal form.
546 *
547  CALL dsytrd( uplo, n, a, lda, work( indd ), work( inde ),
548  $ work( indtau ), work( indwk ), llwork, iinfo )
549 *
550 * If all eigenvalues are desired
551 * then call DSTERF or DSTEMR and DORMTR.
552 *
553  IF( ( alleig .OR. ( indeig .AND. il.EQ.1 .AND. iu.EQ.n ) ) .AND.
554  $ ieeeok.EQ.1 ) THEN
555  IF( .NOT.wantz ) THEN
556  CALL dcopy( n, work( indd ), 1, w, 1 )
557  CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
558  CALL dsterf( n, w, work( indee ), info )
559  ELSE
560  CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
561  CALL dcopy( n, work( indd ), 1, work( inddd ), 1 )
562 *
563  IF (abstol .LE. two*n*eps) THEN
564  tryrac = .true.
565  ELSE
566  tryrac = .false.
567  END IF
568  CALL dstemr( jobz, 'A', n, work( inddd ), work( indee ),
569  $ vl, vu, il, iu, m, w, z, ldz, n, isuppz,
570  $ tryrac, work( indwk ), lwork, iwork, liwork,
571  $ info )
572 *
573 *
574 *
575 * Apply orthogonal matrix used in reduction to tridiagonal
576 * form to eigenvectors returned by DSTEIN.
577 *
578  IF( wantz .AND. info.EQ.0 ) THEN
579  indwkn = inde
580  llwrkn = lwork - indwkn + 1
581  CALL dormtr( 'L', uplo, 'N', n, m, a, lda,
582  $ work( indtau ), z, ldz, work( indwkn ),
583  $ llwrkn, iinfo )
584  END IF
585  END IF
586 *
587 *
588  IF( info.EQ.0 ) THEN
589 * Everything worked. Skip DSTEBZ/DSTEIN. IWORK(:) are
590 * undefined.
591  m = n
592  go to 30
593  END IF
594  info = 0
595  END IF
596 *
597 * Otherwise, call DSTEBZ and, if eigenvectors are desired, DSTEIN.
598 * Also call DSTEBZ and DSTEIN if DSTEMR fails.
599 *
600  IF( wantz ) THEN
601  order = 'B'
602  ELSE
603  order = 'E'
604  END IF
605 
606  CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
607  $ work( indd ), work( inde ), m, nsplit, w,
608  $ iwork( indibl ), iwork( indisp ), work( indwk ),
609  $ iwork( indiwo ), info )
610 *
611  IF( wantz ) THEN
612  CALL dstein( n, work( indd ), work( inde ), m, w,
613  $ iwork( indibl ), iwork( indisp ), z, ldz,
614  $ work( indwk ), iwork( indiwo ), iwork( indifl ),
615  $ info )
616 *
617 * Apply orthogonal matrix used in reduction to tridiagonal
618 * form to eigenvectors returned by DSTEIN.
619 *
620  indwkn = inde
621  llwrkn = lwork - indwkn + 1
622  CALL dormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
623  $ ldz, work( indwkn ), llwrkn, iinfo )
624  END IF
625 *
626 * If matrix was scaled, then rescale eigenvalues appropriately.
627 *
628 * Jump here if DSTEMR/DSTEIN succeeded.
629  30 continue
630  IF( iscale.EQ.1 ) THEN
631  IF( info.EQ.0 ) THEN
632  imax = m
633  ELSE
634  imax = info - 1
635  END IF
636  CALL dscal( imax, one / sigma, w, 1 )
637  END IF
638 *
639 * If eigenvalues are not in order, then sort them, along with
640 * eigenvectors. Note: We do not sort the IFAIL portion of IWORK.
641 * It may not be initialized (if DSTEMR/DSTEIN succeeded), and we do
642 * not return this detailed information to the user.
643 *
644  IF( wantz ) THEN
645  DO 50 j = 1, m - 1
646  i = 0
647  tmp1 = w( j )
648  DO 40 jj = j + 1, m
649  IF( w( jj ).LT.tmp1 ) THEN
650  i = jj
651  tmp1 = w( jj )
652  END IF
653  40 continue
654 *
655  IF( i.NE.0 ) THEN
656  w( i ) = w( j )
657  w( j ) = tmp1
658  CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
659  END IF
660  50 continue
661  END IF
662 *
663 * Set WORK(1) to optimal workspace size.
664 *
665  work( 1 ) = lwkopt
666  iwork( 1 ) = liwmin
667 *
668  return
669 *
670 * End of DSYEVR
671 *
672  END