LAPACK 3.12.0
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*> \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 DSTEIN 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 ). This is an output of DSTEMR (tridiagonal
260*> matrix). The support of the eigenvectors of A is typically
261*> 1:N because of the orthogonal transformations applied by DORMTR.
262*> Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1
263*> \endverbatim
264*>
265*> \param[out] WORK
266*> \verbatim
267*> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
268*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
269*> \endverbatim
270*>
271*> \param[in] LWORK
272*> \verbatim
273*> LWORK is INTEGER
274*> The dimension of the array WORK. LWORK >= max(1,26*N).
275*> For optimal efficiency, LWORK >= (NB+6)*N,
276*> where NB is the max of the blocksize for DSYTRD and DORMTR
277*> returned by ILAENV.
278*>
279*> If LWORK = -1, then a workspace query is assumed; the routine
280*> only calculates the optimal size of the WORK array, returns
281*> this value as the first entry of the WORK array, and no error
282*> message related to LWORK is issued by XERBLA.
283*> \endverbatim
284*>
285*> \param[out] IWORK
286*> \verbatim
287*> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
288*> On exit, if INFO = 0, IWORK(1) returns the optimal LWORK.
289*> \endverbatim
290*>
291*> \param[in] LIWORK
292*> \verbatim
293*> LIWORK is INTEGER
294*> The dimension of the array IWORK. LIWORK >= max(1,10*N).
295*>
296*> If LIWORK = -1, then a workspace query is assumed; the
297*> routine only calculates the optimal size of the IWORK array,
298*> returns this value as the first entry of the IWORK array, and
299*> no error message related to LIWORK is issued by XERBLA.
300*> \endverbatim
301*>
302*> \param[out] INFO
303*> \verbatim
304*> INFO is INTEGER
305*> = 0: successful exit
306*> < 0: if INFO = -i, the i-th argument had an illegal value
307*> > 0: Internal error
308*> \endverbatim
309*
310* Authors:
311* ========
312*
313*> \author Univ. of Tennessee
314*> \author Univ. of California Berkeley
315*> \author Univ. of Colorado Denver
316*> \author NAG Ltd.
317*
318*> \ingroup heevr
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 --
336* -- LAPACK is a software package provided by Univ. of Tennessee, --
337* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
338*
339* .. Scalar Arguments ..
340 CHARACTER JOBZ, RANGE, UPLO
341 INTEGER IL, INFO, IU, LDA, LDZ, LIWORK, LWORK, M, N
342 DOUBLE PRECISION ABSTOL, VL, VU
343* ..
344* .. Array Arguments ..
345 INTEGER ISUPPZ( * ), IWORK( * )
346 DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
347* ..
348*
349* =====================================================================
350*
351* .. Parameters ..
352 DOUBLE PRECISION ZERO, ONE, TWO
353 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0, two = 2.0d+0 )
354* ..
355* .. Local Scalars ..
356 LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, VALEIG, WANTZ,
357 $ TRYRAC
358 CHARACTER ORDER
359 INTEGER I, IEEEOK, IINFO, IMAX, INDD, INDDD, INDE,
360 $ indee, indibl, indifl, indisp, indiwo, indtau,
361 $ indwk, indwkn, iscale, j, jj, liwmin,
362 $ llwork, llwrkn, lwkopt, lwmin, nb, nsplit
363 DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
364 $ SIGMA, SMLNUM, TMP1, VLL, VUU
365* ..
366* .. External Functions ..
367 LOGICAL LSAME
368 INTEGER ILAENV
369 DOUBLE PRECISION DLAMCH, DLANSY
370 EXTERNAL lsame, ilaenv, dlamch, dlansy
371* ..
372* .. External Subroutines ..
373 EXTERNAL dcopy, dormtr, dscal, dstebz, dstemr, dstein,
375* ..
376* .. Intrinsic Functions ..
377 INTRINSIC max, min, sqrt
378* ..
379* .. Executable Statements ..
380*
381* Test the input parameters.
382*
383 ieeeok = ilaenv( 10, 'DSYEVR', 'N', 1, 2, 3, 4 )
384*
385 lower = lsame( uplo, 'L' )
386 wantz = lsame( jobz, 'V' )
387 alleig = lsame( range, 'A' )
388 valeig = lsame( range, 'V' )
389 indeig = lsame( range, 'I' )
390*
391 lquery = ( ( lwork.EQ.-1 ) .OR. ( liwork.EQ.-1 ) )
392*
393 lwmin = max( 1, 26*n )
394 liwmin = max( 1, 10*n )
395*
396 info = 0
397 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
398 info = -1
399 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
400 info = -2
401 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
402 info = -3
403 ELSE IF( n.LT.0 ) THEN
404 info = -4
405 ELSE IF( lda.LT.max( 1, n ) ) THEN
406 info = -6
407 ELSE
408 IF( valeig ) THEN
409 IF( n.GT.0 .AND. vu.LE.vl )
410 $ info = -8
411 ELSE IF( indeig ) THEN
412 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
413 info = -9
414 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
415 info = -10
416 END IF
417 END IF
418 END IF
419 IF( info.EQ.0 ) THEN
420 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
421 info = -15
422 ELSE IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
423 info = -18
424 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
425 info = -20
426 END IF
427 END IF
428*
429 IF( info.EQ.0 ) THEN
430 nb = ilaenv( 1, 'DSYTRD', uplo, n, -1, -1, -1 )
431 nb = max( nb, ilaenv( 1, 'DORMTR', uplo, n, -1, -1, -1 ) )
432 lwkopt = max( ( nb+1 )*n, lwmin )
433 work( 1 ) = lwkopt
434 iwork( 1 ) = liwmin
435 END IF
436*
437 IF( info.NE.0 ) THEN
438 CALL xerbla( 'DSYEVR', -info )
439 RETURN
440 ELSE IF( lquery ) THEN
441 RETURN
442 END IF
443*
444* Quick return if possible
445*
446 m = 0
447 IF( n.EQ.0 ) THEN
448 work( 1 ) = 1
449 RETURN
450 END IF
451*
452 IF( n.EQ.1 ) THEN
453 work( 1 ) = 7
454 IF( alleig .OR. indeig ) THEN
455 m = 1
456 w( 1 ) = a( 1, 1 )
457 ELSE
458 IF( vl.LT.a( 1, 1 ) .AND. vu.GE.a( 1, 1 ) ) THEN
459 m = 1
460 w( 1 ) = a( 1, 1 )
461 END IF
462 END IF
463 IF( wantz ) THEN
464 z( 1, 1 ) = one
465 isuppz( 1 ) = 1
466 isuppz( 2 ) = 1
467 END IF
468 RETURN
469 END IF
470*
471* Get machine constants.
472*
473 safmin = dlamch( 'Safe minimum' )
474 eps = dlamch( 'Precision' )
475 smlnum = safmin / eps
476 bignum = one / smlnum
477 rmin = sqrt( smlnum )
478 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
479*
480* Scale matrix to allowable range, if necessary.
481*
482 iscale = 0
483 abstll = abstol
484 IF (valeig) THEN
485 vll = vl
486 vuu = vu
487 END IF
488 anrm = dlansy( 'M', uplo, n, a, lda, work )
489 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
490 iscale = 1
491 sigma = rmin / anrm
492 ELSE IF( anrm.GT.rmax ) THEN
493 iscale = 1
494 sigma = rmax / anrm
495 END IF
496 IF( iscale.EQ.1 ) THEN
497 IF( lower ) THEN
498 DO 10 j = 1, n
499 CALL dscal( n-j+1, sigma, a( j, j ), 1 )
500 10 CONTINUE
501 ELSE
502 DO 20 j = 1, n
503 CALL dscal( j, sigma, a( 1, j ), 1 )
504 20 CONTINUE
505 END IF
506 IF( abstol.GT.0 )
507 $ abstll = abstol*sigma
508 IF( valeig ) THEN
509 vll = vl*sigma
510 vuu = vu*sigma
511 END IF
512 END IF
513
514* Initialize indices into workspaces. Note: The IWORK indices are
515* used only if DSTERF or DSTEMR fail.
516
517* WORK(INDTAU:INDTAU+N-1) stores the scalar factors of the
518* elementary reflectors used in DSYTRD.
519 indtau = 1
520* WORK(INDD:INDD+N-1) stores the tridiagonal's diagonal entries.
521 indd = indtau + n
522* WORK(INDE:INDE+N-1) stores the off-diagonal entries of the
523* tridiagonal matrix from DSYTRD.
524 inde = indd + n
525* WORK(INDDD:INDDD+N-1) is a copy of the diagonal entries over
526* -written by DSTEMR (the DSTERF path copies the diagonal to W).
527 inddd = inde + n
528* WORK(INDEE:INDEE+N-1) is a copy of the off-diagonal entries over
529* -written while computing the eigenvalues in DSTERF and DSTEMR.
530 indee = inddd + n
531* INDWK is the starting offset of the left-over workspace, and
532* LLWORK is the remaining workspace size.
533 indwk = indee + n
534 llwork = lwork - indwk + 1
535
536* IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in DSTEBZ and
537* stores the block indices of each of the M<=N eigenvalues.
538 indibl = 1
539* IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in DSTEBZ and
540* stores the starting and finishing indices of each block.
541 indisp = indibl + n
542* IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
543* that corresponding to eigenvectors that fail to converge in
544* DSTEIN. This information is discarded; if any fail, the driver
545* returns INFO > 0.
546 indifl = indisp + n
547* INDIWO is the offset of the remaining integer workspace.
548 indiwo = indifl + n
549
550*
551* Call DSYTRD to reduce symmetric matrix to tridiagonal form.
552*
553 CALL dsytrd( uplo, n, a, lda, work( indd ), work( inde ),
554 $ work( indtau ), work( indwk ), llwork, iinfo )
555*
556* If all eigenvalues are desired
557* then call DSTERF or DSTEMR and DORMTR.
558*
559 IF( ( alleig .OR. ( indeig .AND. il.EQ.1 .AND. iu.EQ.n ) ) .AND.
560 $ ieeeok.EQ.1 ) THEN
561 IF( .NOT.wantz ) THEN
562 CALL dcopy( n, work( indd ), 1, w, 1 )
563 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
564 CALL dsterf( n, w, work( indee ), info )
565 ELSE
566 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
567 CALL dcopy( n, work( indd ), 1, work( inddd ), 1 )
568*
569 IF (abstol .LE. two*n*eps) THEN
570 tryrac = .true.
571 ELSE
572 tryrac = .false.
573 END IF
574 CALL dstemr( jobz, 'A', n, work( inddd ), work( indee ),
575 $ vl, vu, il, iu, m, w, z, ldz, n, isuppz,
576 $ tryrac, work( indwk ), lwork, iwork, liwork,
577 $ info )
578*
579*
580*
581* Apply orthogonal matrix used in reduction to tridiagonal
582* form to eigenvectors returned by DSTEMR.
583*
584 IF( wantz .AND. info.EQ.0 ) THEN
585 indwkn = inde
586 llwrkn = lwork - indwkn + 1
587 CALL dormtr( 'L', uplo, 'N', n, m, a, lda,
588 $ work( indtau ), z, ldz, work( indwkn ),
589 $ llwrkn, iinfo )
590 END IF
591 END IF
592*
593*
594 IF( info.EQ.0 ) THEN
595* Everything worked. Skip DSTEBZ/DSTEIN. IWORK(:) are
596* undefined.
597 m = n
598 GO TO 30
599 END IF
600 info = 0
601 END IF
602*
603* Otherwise, call DSTEBZ and, if eigenvectors are desired, DSTEIN.
604* Also call DSTEBZ and DSTEIN if DSTEMR fails.
605*
606 IF( wantz ) THEN
607 order = 'B'
608 ELSE
609 order = 'E'
610 END IF
611
612 CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
613 $ work( indd ), work( inde ), m, nsplit, w,
614 $ iwork( indibl ), iwork( indisp ), work( indwk ),
615 $ iwork( indiwo ), info )
616*
617 IF( wantz ) THEN
618 CALL dstein( n, work( indd ), work( inde ), m, w,
619 $ iwork( indibl ), iwork( indisp ), z, ldz,
620 $ work( indwk ), iwork( indiwo ), iwork( indifl ),
621 $ info )
622*
623* Apply orthogonal matrix used in reduction to tridiagonal
624* form to eigenvectors returned by DSTEIN.
625*
626 indwkn = inde
627 llwrkn = lwork - indwkn + 1
628 CALL dormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
629 $ ldz, work( indwkn ), llwrkn, iinfo )
630 END IF
631*
632* If matrix was scaled, then rescale eigenvalues appropriately.
633*
634* Jump here if DSTEMR/DSTEIN succeeded.
635 30 CONTINUE
636 IF( iscale.EQ.1 ) THEN
637 IF( info.EQ.0 ) THEN
638 imax = m
639 ELSE
640 imax = info - 1
641 END IF
642 CALL dscal( imax, one / sigma, w, 1 )
643 END IF
644*
645* If eigenvalues are not in order, then sort them, along with
646* eigenvectors. Note: We do not sort the IFAIL portion of IWORK.
647* It may not be initialized (if DSTEMR/DSTEIN succeeded), and we do
648* not return this detailed information to the user.
649*
650 IF( wantz ) THEN
651 DO 50 j = 1, m - 1
652 i = 0
653 tmp1 = w( j )
654 DO 40 jj = j + 1, m
655 IF( w( jj ).LT.tmp1 ) THEN
656 i = jj
657 tmp1 = w( jj )
658 END IF
659 40 CONTINUE
660*
661 IF( i.NE.0 ) THEN
662 w( i ) = w( j )
663 w( j ) = tmp1
664 CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
665 END IF
666 50 CONTINUE
667 END IF
668*
669* Set WORK(1) to optimal workspace size.
670*
671 work( 1 ) = lwkopt
672 iwork( 1 ) = liwmin
673*
674 RETURN
675*
676* End of DSYEVR
677*
678 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:334
subroutine dsytrd(uplo, n, a, lda, d, e, tau, work, lwork, info)
DSYTRD
Definition dsytrd.f:192
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:273
subroutine dstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
DSTEIN
Definition dstein.f:174
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:322
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
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:171