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