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