LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
sstevr.f
Go to the documentation of this file.
1*> \brief <b> SSTEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SSTEVR + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sstevr.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sstevr.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sstevr.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SSTEVR( JOBZ, RANGE, N, D, E, VL, VU, IL, IU, ABSTOL,
20* M, W, Z, LDZ, ISUPPZ, WORK, LWORK, IWORK,
21* LIWORK, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER JOBZ, RANGE
25* INTEGER IL, INFO, IU, LDZ, LIWORK, LWORK, M, N
26* REAL ABSTOL, VL, VU
27* ..
28* .. Array Arguments ..
29* INTEGER ISUPPZ( * ), IWORK( * )
30* REAL D( * ), E( * ), W( * ), WORK( * ), Z( LDZ, * )
31* ..
32*
33*
34*> \par Purpose:
35* =============
36*>
37*> \verbatim
38*>
39*> SSTEVR computes selected eigenvalues and, optionally, eigenvectors
40*> of a real symmetric tridiagonal matrix T. Eigenvalues and
41*> eigenvectors can be selected by specifying either a range of values
42*> or a range of indices for the desired eigenvalues.
43*>
44*> Whenever possible, SSTEVR calls SSTEMR to compute the
45*> eigenspectrum using Relatively Robust Representations. SSTEMR
46*> computes eigenvalues by the dqds algorithm, while orthogonal
47*> eigenvectors are computed from various "good" L D L^T representations
48*> (also known as Relatively Robust Representations). Gram-Schmidt
49*> orthogonalization is avoided as far as possible. More specifically,
50*> the various steps of the algorithm are as follows. For the i-th
51*> unreduced block of T,
52*> (a) Compute T - sigma_i = L_i D_i L_i^T, such that L_i D_i L_i^T
53*> is a relatively robust representation,
54*> (b) Compute the eigenvalues, lambda_j, of L_i D_i L_i^T to high
55*> relative accuracy by the dqds algorithm,
56*> (c) If there is a cluster of close eigenvalues, "choose" sigma_i
57*> close to the cluster, and go to step (a),
58*> (d) Given the approximate eigenvalue lambda_j of L_i D_i L_i^T,
59*> compute the corresponding eigenvector by forming a
60*> rank-revealing twisted factorization.
61*> The desired accuracy of the output can be specified by the input
62*> parameter ABSTOL.
63*>
64*> For more details, see "A new O(n^2) algorithm for the symmetric
65*> tridiagonal eigenvalue/eigenvector problem", by Inderjit Dhillon,
66*> Computer Science Division Technical Report No. UCB//CSD-97-971,
67*> UC Berkeley, May 1997.
68*>
69*>
70*> Note 1 : SSTEVR calls SSTEMR when the full spectrum is requested
71*> on machines which conform to the ieee-754 floating point standard.
72*> SSTEVR calls SSTEBZ and SSTEIN on non-ieee machines and
73*> when partial spectrum requests are made.
74*>
75*> Normal execution of SSTEMR may create NaNs and infinities and
76*> hence may abort due to a floating point exception in environments
77*> which do not handle NaNs and infinities in the ieee standard default
78*> manner.
79*> \endverbatim
80*
81* Arguments:
82* ==========
83*
84*> \param[in] JOBZ
85*> \verbatim
86*> JOBZ is CHARACTER*1
87*> = 'N': Compute eigenvalues only;
88*> = 'V': Compute eigenvalues and eigenvectors.
89*> \endverbatim
90*>
91*> \param[in] RANGE
92*> \verbatim
93*> RANGE is CHARACTER*1
94*> = 'A': all eigenvalues will be found.
95*> = 'V': all eigenvalues in the half-open interval (VL,VU]
96*> will be found.
97*> = 'I': the IL-th through IU-th eigenvalues will be found.
98*> For RANGE = 'V' or 'I' and IU - IL < N - 1, SSTEBZ and
99*> SSTEIN are called
100*> \endverbatim
101*>
102*> \param[in] N
103*> \verbatim
104*> N is INTEGER
105*> The order of the matrix. N >= 0.
106*> \endverbatim
107*>
108*> \param[in,out] D
109*> \verbatim
110*> D is REAL array, dimension (N)
111*> On entry, the n diagonal elements of the tridiagonal matrix
112*> A.
113*> On exit, D may be multiplied by a constant factor chosen
114*> to avoid over/underflow in computing the eigenvalues.
115*> \endverbatim
116*>
117*> \param[in,out] E
118*> \verbatim
119*> E is REAL array, dimension (max(1,N-1))
120*> On entry, the (n-1) subdiagonal elements of the tridiagonal
121*> matrix A in elements 1 to N-1 of E.
122*> On exit, E may be multiplied by a constant factor chosen
123*> to avoid over/underflow in computing the eigenvalues.
124*> \endverbatim
125*>
126*> \param[in] VL
127*> \verbatim
128*> VL is REAL
129*> If RANGE='V', the lower bound of the interval to
130*> be searched for eigenvalues. VL < VU.
131*> Not referenced if RANGE = 'A' or 'I'.
132*> \endverbatim
133*>
134*> \param[in] VU
135*> \verbatim
136*> VU is REAL
137*> If RANGE='V', the upper bound of the interval to
138*> be searched for eigenvalues. VL < VU.
139*> Not referenced if RANGE = 'A' or 'I'.
140*> \endverbatim
141*>
142*> \param[in] IL
143*> \verbatim
144*> IL is INTEGER
145*> If RANGE='I', the index of the
146*> smallest eigenvalue to be returned.
147*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
148*> Not referenced if RANGE = 'A' or 'V'.
149*> \endverbatim
150*>
151*> \param[in] IU
152*> \verbatim
153*> IU is INTEGER
154*> If RANGE='I', the index of the
155*> largest eigenvalue to be returned.
156*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
157*> Not referenced if RANGE = 'A' or 'V'.
158*> \endverbatim
159*>
160*> \param[in] ABSTOL
161*> \verbatim
162*> ABSTOL is REAL
163*> The absolute error tolerance for the eigenvalues.
164*> An approximate eigenvalue is accepted as converged
165*> when it is determined to lie in an interval [a,b]
166*> of width less than or equal to
167*>
168*> ABSTOL + EPS * max( |a|,|b| ) ,
169*>
170*> where EPS is the machine precision. If ABSTOL is less than
171*> or equal to zero, then EPS*|T| will be used in its place,
172*> where |T| is the 1-norm of the tridiagonal matrix obtained
173*> by reducing A to tridiagonal form.
174*>
175*> See "Computing Small Singular Values of Bidiagonal Matrices
176*> with Guaranteed High Relative Accuracy," by Demmel and
177*> Kahan, LAPACK Working Note #3.
178*>
179*> If high relative accuracy is important, set ABSTOL to
180*> SLAMCH( 'Safe minimum' ). Doing so will guarantee that
181*> eigenvalues are computed to high relative accuracy when
182*> possible in future releases. The current code does not
183*> make any guarantees about high relative accuracy, but
184*> future releases will. See J. Barlow and J. Demmel,
185*> "Computing Accurate Eigensystems of Scaled Diagonally
186*> Dominant Matrices", LAPACK Working Note #7, for a discussion
187*> of which matrices define their eigenvalues to high relative
188*> accuracy.
189*> \endverbatim
190*>
191*> \param[out] M
192*> \verbatim
193*> M is INTEGER
194*> The total number of eigenvalues found. 0 <= M <= N.
195*> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
196*> \endverbatim
197*>
198*> \param[out] W
199*> \verbatim
200*> W is REAL array, dimension (N)
201*> The first M elements contain the selected eigenvalues in
202*> ascending order.
203*> \endverbatim
204*>
205*> \param[out] Z
206*> \verbatim
207*> Z is REAL array, dimension (LDZ, max(1,M) )
208*> If JOBZ = 'V', then if INFO = 0, the first M columns of Z
209*> contain the orthonormal eigenvectors of the matrix A
210*> corresponding to the selected eigenvalues, with the i-th
211*> column of Z holding the eigenvector associated with W(i).
212*> Note: the user must ensure that at least max(1,M) columns are
213*> supplied in the array Z; if RANGE = 'V', the exact value of M
214*> is not known in advance and an upper bound must be used.
215*> \endverbatim
216*>
217*> \param[in] LDZ
218*> \verbatim
219*> LDZ is INTEGER
220*> The leading dimension of the array Z. LDZ >= 1, and if
221*> JOBZ = 'V', LDZ >= max(1,N).
222*> \endverbatim
223*>
224*> \param[out] ISUPPZ
225*> \verbatim
226*> ISUPPZ is INTEGER array, dimension ( 2*max(1,M) )
227*> The support of the eigenvectors in Z, i.e., the indices
228*> indicating the nonzero elements in Z. The i-th eigenvector
229*> is nonzero only in elements ISUPPZ( 2*i-1 ) through
230*> ISUPPZ( 2*i ).
231*> Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1
232*> \endverbatim
233*>
234*> \param[out] WORK
235*> \verbatim
236*> WORK is REAL array, dimension (MAX(1,LWORK))
237*> On exit, if INFO = 0, WORK(1) returns the optimal (and
238*> minimal) LWORK.
239*> \endverbatim
240*>
241*> \param[in] LWORK
242*> \verbatim
243*> LWORK is INTEGER
244*> The dimension of the array WORK. LWORK >= 20*N.
245*>
246*> If LWORK = -1, then a workspace query is assumed; the routine
247*> only calculates the optimal sizes of the WORK and IWORK
248*> arrays, returns these values as the first entries of the WORK
249*> and IWORK arrays, and no error message related to LWORK or
250*> LIWORK is issued by XERBLA.
251*> \endverbatim
252*>
253*> \param[out] IWORK
254*> \verbatim
255*> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
256*> On exit, if INFO = 0, IWORK(1) returns the optimal (and
257*> minimal) LIWORK.
258*> \endverbatim
259*>
260*> \param[in] LIWORK
261*> \verbatim
262*> LIWORK is INTEGER
263*> The dimension of the array IWORK. LIWORK >= 10*N.
264*>
265*> If LIWORK = -1, then a workspace query is assumed; the
266*> routine only calculates the optimal sizes of the WORK and
267*> IWORK arrays, returns these values as the first entries of
268*> the WORK and IWORK arrays, and no error message related to
269*> LWORK or LIWORK is issued by XERBLA.
270*> \endverbatim
271*>
272*> \param[out] INFO
273*> \verbatim
274*> INFO is INTEGER
275*> = 0: successful exit
276*> < 0: if INFO = -i, the i-th argument had an illegal value
277*> > 0: Internal error
278*> \endverbatim
279*
280* Authors:
281* ========
282*
283*> \author Univ. of Tennessee
284*> \author Univ. of California Berkeley
285*> \author Univ. of Colorado Denver
286*> \author NAG Ltd.
287*
288*> \ingroup stevr
289*
290*> \par Contributors:
291* ==================
292*>
293*> Inderjit Dhillon, IBM Almaden, USA \n
294*> Osni Marques, LBNL/NERSC, USA \n
295*> Ken Stanley, Computer Science Division, University of
296*> California at Berkeley, USA \n
297*> Jason Riedy, Computer Science Division, University of
298*> California at Berkeley, USA \n
299*>
300* =====================================================================
301 SUBROUTINE sstevr( JOBZ, RANGE, N, D, E, VL, VU, IL, IU,
302 $ ABSTOL,
303 $ M, W, Z, LDZ, ISUPPZ, WORK, LWORK, IWORK,
304 $ LIWORK, INFO )
305*
306* -- LAPACK driver routine --
307* -- LAPACK is a software package provided by Univ. of Tennessee, --
308* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
309*
310* .. Scalar Arguments ..
311 CHARACTER JOBZ, RANGE
312 INTEGER IL, INFO, IU, LDZ, LIWORK, LWORK, M, N
313 REAL ABSTOL, VL, VU
314* ..
315* .. Array Arguments ..
316 INTEGER ISUPPZ( * ), IWORK( * )
317 REAL D( * ), E( * ), W( * ), WORK( * ), Z( LDZ, * )
318* ..
319*
320* =====================================================================
321*
322* .. Parameters ..
323 REAL ZERO, ONE, TWO
324 PARAMETER ( ZERO = 0.0e+0, one = 1.0e+0, two = 2.0e+0 )
325* ..
326* .. Local Scalars ..
327 LOGICAL ALLEIG, INDEIG, TEST, LQUERY, VALEIG, WANTZ,
328 $ TRYRAC
329 CHARACTER ORDER
330 INTEGER I, IEEEOK, IMAX, INDIBL, INDIFL, INDISP,
331 $ INDIWO, ISCALE, J, JJ, LIWMIN, LWMIN, NSPLIT
332 REAL BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA, SMLNUM,
333 $ tmp1, tnrm, vll, vuu
334* ..
335* .. External Functions ..
336 LOGICAL LSAME
337 INTEGER ILAENV
338 REAL SLAMCH, SLANST, SROUNDUP_LWORK
339 EXTERNAL lsame, ilaenv, slamch, slanst,
340 $ sroundup_lwork
341* ..
342* .. External Subroutines ..
343 EXTERNAL scopy, sscal, sstebz, sstemr, sstein,
344 $ ssterf,
345 $ sswap, xerbla
346* ..
347* .. Intrinsic Functions ..
348 INTRINSIC max, min, sqrt
349* ..
350* .. Executable Statements ..
351*
352*
353* Test the input parameters.
354*
355 ieeeok = ilaenv( 10, 'SSTEVR', 'N', 1, 2, 3, 4 )
356*
357 wantz = lsame( jobz, 'V' )
358 alleig = lsame( range, 'A' )
359 valeig = lsame( range, 'V' )
360 indeig = lsame( range, 'I' )
361*
362 lquery = ( ( lwork.EQ.-1 ) .OR. ( liwork.EQ.-1 ) )
363 lwmin = max( 1, 20*n )
364 liwmin = max(1, 10*n )
365*
366*
367 info = 0
368 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
369 info = -1
370 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
371 info = -2
372 ELSE IF( n.LT.0 ) THEN
373 info = -3
374 ELSE
375 IF( valeig ) THEN
376 IF( n.GT.0 .AND. vu.LE.vl )
377 $ info = -7
378 ELSE IF( indeig ) THEN
379 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
380 info = -8
381 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
382 info = -9
383 END IF
384 END IF
385 END IF
386 IF( info.EQ.0 ) THEN
387 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
388 info = -14
389 END IF
390 END IF
391*
392 IF( info.EQ.0 ) THEN
393 work( 1 ) = sroundup_lwork(lwmin)
394 iwork( 1 ) = liwmin
395*
396 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
397 info = -17
398 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
399 info = -19
400 END IF
401 END IF
402*
403 IF( info.NE.0 ) THEN
404 CALL xerbla( 'SSTEVR', -info )
405 RETURN
406 ELSE IF( lquery ) THEN
407 RETURN
408 END IF
409*
410* Quick return if possible
411*
412 m = 0
413 IF( n.EQ.0 )
414 $ RETURN
415*
416 IF( n.EQ.1 ) THEN
417 IF( alleig .OR. indeig ) THEN
418 m = 1
419 w( 1 ) = d( 1 )
420 ELSE
421 IF( vl.LT.d( 1 ) .AND. vu.GE.d( 1 ) ) THEN
422 m = 1
423 w( 1 ) = d( 1 )
424 END IF
425 END IF
426 IF( wantz )
427 $ z( 1, 1 ) = one
428 RETURN
429 END IF
430*
431* Get machine constants.
432*
433 safmin = slamch( 'Safe minimum' )
434 eps = slamch( 'Precision' )
435 smlnum = safmin / eps
436 bignum = one / smlnum
437 rmin = sqrt( smlnum )
438 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
439*
440*
441* Scale matrix to allowable range, if necessary.
442*
443 iscale = 0
444 IF( valeig ) THEN
445 vll = vl
446 vuu = vu
447 END IF
448*
449 tnrm = slanst( 'M', n, d, e )
450 IF( tnrm.GT.zero .AND. tnrm.LT.rmin ) THEN
451 iscale = 1
452 sigma = rmin / tnrm
453 ELSE IF( tnrm.GT.rmax ) THEN
454 iscale = 1
455 sigma = rmax / tnrm
456 END IF
457 IF( iscale.EQ.1 ) THEN
458 CALL sscal( n, sigma, d, 1 )
459 CALL sscal( n-1, sigma, e( 1 ), 1 )
460 IF( valeig ) THEN
461 vll = vl*sigma
462 vuu = vu*sigma
463 END IF
464 END IF
465
466* Initialize indices into workspaces. Note: These indices are used only
467* if SSTERF or SSTEMR fail.
468
469* IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in SSTEBZ and
470* stores the block indices of each of the M<=N eigenvalues.
471 indibl = 1
472* IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in SSTEBZ and
473* stores the starting and finishing indices of each block.
474 indisp = indibl + n
475* IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
476* that corresponding to eigenvectors that fail to converge in
477* SSTEIN. This information is discarded; if any fail, the driver
478* returns INFO > 0.
479 indifl = indisp + n
480* INDIWO is the offset of the remaining integer workspace.
481 indiwo = indisp + n
482*
483* If all eigenvalues are desired, then
484* call SSTERF or SSTEMR. If this fails for some eigenvalue, then
485* try SSTEBZ.
486*
487*
488 test = .false.
489 IF( indeig ) THEN
490 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
491 test = .true.
492 END IF
493 END IF
494 IF( ( alleig .OR. test ) .AND. ieeeok.EQ.1 ) THEN
495 CALL scopy( n-1, e( 1 ), 1, work( 1 ), 1 )
496 IF( .NOT.wantz ) THEN
497 CALL scopy( n, d, 1, w, 1 )
498 CALL ssterf( n, w, work, info )
499 ELSE
500 CALL scopy( n, d, 1, work( n+1 ), 1 )
501 IF (abstol .LE. two*real( n )*eps) THEN
502 tryrac = .true.
503 ELSE
504 tryrac = .false.
505 END IF
506 CALL sstemr( jobz, 'A', n, work( n+1 ), work, vl, vu, il,
507 $ iu, m, w, z, ldz, n, isuppz, tryrac,
508 $ work( 2*n+1 ), lwork-2*n, iwork, liwork, info )
509*
510 END IF
511 IF( info.EQ.0 ) THEN
512 m = n
513 GO TO 10
514 END IF
515 info = 0
516 END IF
517*
518* Otherwise, call SSTEBZ and, if eigenvectors are desired, SSTEIN.
519*
520 IF( wantz ) THEN
521 order = 'B'
522 ELSE
523 order = 'E'
524 END IF
525
526 CALL sstebz( range, order, n, vll, vuu, il, iu, abstol, d, e,
527 $ m,
528 $ nsplit, w, iwork( indibl ), iwork( indisp ), work,
529 $ iwork( indiwo ), info )
530*
531 IF( wantz ) THEN
532 CALL sstein( n, d, e, m, w, iwork( indibl ),
533 $ iwork( indisp ),
534 $ z, ldz, work, iwork( indiwo ), iwork( indifl ),
535 $ info )
536 END IF
537*
538* If matrix was scaled, then rescale eigenvalues appropriately.
539*
540 10 CONTINUE
541 IF( iscale.EQ.1 ) THEN
542 IF( info.EQ.0 ) THEN
543 imax = m
544 ELSE
545 imax = info - 1
546 END IF
547 CALL sscal( imax, one / sigma, w, 1 )
548 END IF
549*
550* If eigenvalues are not in order, then sort them, along with
551* eigenvectors.
552*
553 IF( wantz ) THEN
554 DO 30 j = 1, m - 1
555 i = 0
556 tmp1 = w( j )
557 DO 20 jj = j + 1, m
558 IF( w( jj ).LT.tmp1 ) THEN
559 i = jj
560 tmp1 = w( jj )
561 END IF
562 20 CONTINUE
563*
564 IF( i.NE.0 ) THEN
565 w( i ) = w( j )
566 w( j ) = tmp1
567 CALL sswap( n, z( 1, i ), 1, z( 1, j ), 1 )
568 END IF
569 30 CONTINUE
570 END IF
571*
572* Causes problems with tests 19 & 20:
573* IF (wantz .and. INDEIG ) Z( 1,1) = Z(1,1) / 1.002 + .002
574*
575*
576 work( 1 ) = sroundup_lwork(lwmin)
577 iwork( 1 ) = liwmin
578 RETURN
579*
580* End of SSTEVR
581*
582 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine sstebz(range, order, n, vl, vu, il, iu, abstol, d, e, m, nsplit, w, iblock, isplit, work, iwork, info)
SSTEBZ
Definition sstebz.f:272
subroutine sstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
SSTEIN
Definition sstein.f:172
subroutine sstemr(jobz, range, n, d, e, vl, vu, il, iu, m, w, z, ldz, nzc, isuppz, tryrac, work, lwork, iwork, liwork, info)
SSTEMR
Definition sstemr.f:320
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:84
subroutine sstevr(jobz, range, n, d, e, vl, vu, il, iu, abstol, m, w, z, ldz, isuppz, work, lwork, iwork, liwork, info)
SSTEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition sstevr.f:305
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82