LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dstevr.f
Go to the documentation of this file.
1*> \brief <b> DSTEVR 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 DSTEVR + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dstevr.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dstevr.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dstevr.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DSTEVR( 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* DOUBLE PRECISION ABSTOL, VL, VU
27* ..
28* .. Array Arguments ..
29* INTEGER ISUPPZ( * ), IWORK( * )
30* DOUBLE PRECISION D( * ), E( * ), W( * ), WORK( * ), Z( LDZ, * )
31* ..
32*
33*
34*> \par Purpose:
35* =============
36*>
37*> \verbatim
38*>
39*> DSTEVR 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, DSTEVR calls DSTEMR to compute the
45*> eigenspectrum using Relatively Robust Representations. DSTEMR
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 : DSTEVR calls DSTEMR when the full spectrum is requested
71*> on machines which conform to the ieee-754 floating point standard.
72*> DSTEVR calls DSTEBZ and DSTEIN on non-ieee machines and
73*> when partial spectrum requests are made.
74*>
75*> Normal execution of DSTEMR 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, DSTEBZ and
99*> DSTEIN 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION
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 DOUBLE PRECISION
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 DOUBLE PRECISION
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*> DLAMCH( '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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 >= max(1,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 >= max(1,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*>
298* =====================================================================
299 SUBROUTINE dstevr( JOBZ, RANGE, N, D, E, VL, VU, IL, IU,
300 $ ABSTOL,
301 $ M, W, Z, LDZ, ISUPPZ, WORK, LWORK, IWORK,
302 $ LIWORK, INFO )
303*
304* -- LAPACK driver routine --
305* -- LAPACK is a software package provided by Univ. of Tennessee, --
306* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
307*
308* .. Scalar Arguments ..
309 CHARACTER JOBZ, RANGE
310 INTEGER IL, INFO, IU, LDZ, LIWORK, LWORK, M, N
311 DOUBLE PRECISION ABSTOL, VL, VU
312* ..
313* .. Array Arguments ..
314 INTEGER ISUPPZ( * ), IWORK( * )
315 DOUBLE PRECISION D( * ), E( * ), W( * ), WORK( * ), Z( LDZ, * )
316* ..
317*
318* =====================================================================
319*
320* .. Parameters ..
321 DOUBLE PRECISION ZERO, ONE, TWO
322 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0, two = 2.0d+0 )
323* ..
324* .. Local Scalars ..
325 LOGICAL ALLEIG, INDEIG, TEST, LQUERY, VALEIG, WANTZ,
326 $ TRYRAC
327 CHARACTER ORDER
328 INTEGER I, IEEEOK, IMAX, INDIBL, INDIFL, INDISP,
329 $ INDIWO, ISCALE, ITMP1, J, JJ, LIWMIN, LWMIN,
330 $ nsplit
331 DOUBLE PRECISION BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA, SMLNUM,
332 $ TMP1, TNRM, VLL, VUU
333* ..
334* .. External Functions ..
335 LOGICAL LSAME
336 INTEGER ILAENV
337 DOUBLE PRECISION DLAMCH, DLANST
338 EXTERNAL lsame, ilaenv, dlamch, dlanst
339* ..
340* .. External Subroutines ..
341 EXTERNAL dcopy, dscal, dstebz, dstemr, dstein,
342 $ dsterf,
343 $ dswap, xerbla
344* ..
345* .. Intrinsic Functions ..
346 INTRINSIC max, min, sqrt
347* ..
348* .. Executable Statements ..
349*
350*
351* Test the input parameters.
352*
353 ieeeok = ilaenv( 10, 'DSTEVR', 'N', 1, 2, 3, 4 )
354*
355 wantz = lsame( jobz, 'V' )
356 alleig = lsame( range, 'A' )
357 valeig = lsame( range, 'V' )
358 indeig = lsame( range, 'I' )
359*
360 lquery = ( ( lwork.EQ.-1 ) .OR. ( liwork.EQ.-1 ) )
361 lwmin = max( 1, 20*n )
362 liwmin = max( 1, 10*n )
363*
364*
365 info = 0
366 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
367 info = -1
368 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
369 info = -2
370 ELSE IF( n.LT.0 ) THEN
371 info = -3
372 ELSE
373 IF( valeig ) THEN
374 IF( n.GT.0 .AND. vu.LE.vl )
375 $ info = -7
376 ELSE IF( indeig ) THEN
377 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
378 info = -8
379 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
380 info = -9
381 END IF
382 END IF
383 END IF
384 IF( info.EQ.0 ) THEN
385 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
386 info = -14
387 END IF
388 END IF
389*
390 IF( info.EQ.0 ) THEN
391 work( 1 ) = lwmin
392 iwork( 1 ) = liwmin
393*
394 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
395 info = -17
396 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
397 info = -19
398 END IF
399 END IF
400*
401 IF( info.NE.0 ) THEN
402 CALL xerbla( 'DSTEVR', -info )
403 RETURN
404 ELSE IF( lquery ) THEN
405 RETURN
406 END IF
407*
408* Quick return if possible
409*
410 m = 0
411 IF( n.EQ.0 )
412 $ RETURN
413*
414 IF( n.EQ.1 ) THEN
415 IF( alleig .OR. indeig ) THEN
416 m = 1
417 w( 1 ) = d( 1 )
418 ELSE
419 IF( vl.LT.d( 1 ) .AND. vu.GE.d( 1 ) ) THEN
420 m = 1
421 w( 1 ) = d( 1 )
422 END IF
423 END IF
424 IF( wantz )
425 $ z( 1, 1 ) = one
426 RETURN
427 END IF
428*
429* Get machine constants.
430*
431 safmin = dlamch( 'Safe minimum' )
432 eps = dlamch( 'Precision' )
433 smlnum = safmin / eps
434 bignum = one / smlnum
435 rmin = sqrt( smlnum )
436 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
437*
438*
439* Scale matrix to allowable range, if necessary.
440*
441 iscale = 0
442 IF( valeig ) THEN
443 vll = vl
444 vuu = vu
445 END IF
446*
447 tnrm = dlanst( 'M', n, d, e )
448 IF( tnrm.GT.zero .AND. tnrm.LT.rmin ) THEN
449 iscale = 1
450 sigma = rmin / tnrm
451 ELSE IF( tnrm.GT.rmax ) THEN
452 iscale = 1
453 sigma = rmax / tnrm
454 END IF
455 IF( iscale.EQ.1 ) THEN
456 CALL dscal( n, sigma, d, 1 )
457 CALL dscal( n-1, sigma, e( 1 ), 1 )
458 IF( valeig ) THEN
459 vll = vl*sigma
460 vuu = vu*sigma
461 END IF
462 END IF
463
464* Initialize indices into workspaces. Note: These indices are used only
465* if DSTERF or DSTEMR fail.
466
467* IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in DSTEBZ and
468* stores the block indices of each of the M<=N eigenvalues.
469 indibl = 1
470* IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in DSTEBZ and
471* stores the starting and finishing indices of each block.
472 indisp = indibl + n
473* IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors
474* that corresponding to eigenvectors that fail to converge in
475* DSTEIN. This information is discarded; if any fail, the driver
476* returns INFO > 0.
477 indifl = indisp + n
478* INDIWO is the offset of the remaining integer workspace.
479 indiwo = indisp + n
480*
481* If all eigenvalues are desired, then
482* call DSTERF or DSTEMR. If this fails for some eigenvalue, then
483* try DSTEBZ.
484*
485*
486 test = .false.
487 IF( indeig ) THEN
488 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
489 test = .true.
490 END IF
491 END IF
492 IF( ( alleig .OR. test ) .AND. ieeeok.EQ.1 ) THEN
493 CALL dcopy( n-1, e( 1 ), 1, work( 1 ), 1 )
494 IF( .NOT.wantz ) THEN
495 CALL dcopy( n, d, 1, w, 1 )
496 CALL dsterf( n, w, work, info )
497 ELSE
498 CALL dcopy( n, d, 1, work( n+1 ), 1 )
499 IF (abstol .LE. two*n*eps) THEN
500 tryrac = .true.
501 ELSE
502 tryrac = .false.
503 END IF
504 CALL dstemr( jobz, 'A', n, work( n+1 ), work, vl, vu, il,
505 $ iu, m, w, z, ldz, n, isuppz, tryrac,
506 $ work( 2*n+1 ), lwork-2*n, iwork, liwork, info )
507*
508 END IF
509 IF( info.EQ.0 ) THEN
510 m = n
511 GO TO 10
512 END IF
513 info = 0
514 END IF
515*
516* Otherwise, call DSTEBZ and, if eigenvectors are desired, DSTEIN.
517*
518 IF( wantz ) THEN
519 order = 'B'
520 ELSE
521 order = 'E'
522 END IF
523
524 CALL dstebz( range, order, n, vll, vuu, il, iu, abstol, d, e,
525 $ m,
526 $ nsplit, w, iwork( indibl ), iwork( indisp ), work,
527 $ iwork( indiwo ), info )
528*
529 IF( wantz ) THEN
530 CALL dstein( n, d, e, m, w, iwork( indibl ),
531 $ iwork( indisp ),
532 $ z, ldz, work, iwork( indiwo ), iwork( indifl ),
533 $ info )
534 END IF
535*
536* If matrix was scaled, then rescale eigenvalues appropriately.
537*
538 10 CONTINUE
539 IF( iscale.EQ.1 ) THEN
540 IF( info.EQ.0 ) THEN
541 imax = m
542 ELSE
543 imax = info - 1
544 END IF
545 CALL dscal( imax, one / sigma, w, 1 )
546 END IF
547*
548* If eigenvalues are not in order, then sort them, along with
549* eigenvectors.
550*
551 IF( wantz ) THEN
552 DO 30 j = 1, m - 1
553 i = 0
554 tmp1 = w( j )
555 DO 20 jj = j + 1, m
556 IF( w( jj ).LT.tmp1 ) THEN
557 i = jj
558 tmp1 = w( jj )
559 END IF
560 20 CONTINUE
561*
562 IF( i.NE.0 ) THEN
563 itmp1 = iwork( i )
564 w( i ) = w( j )
565 iwork( i ) = iwork( j )
566 w( j ) = tmp1
567 iwork( j ) = itmp1
568 CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
569 END IF
570 30 CONTINUE
571 END IF
572*
573* Causes problems with tests 19 & 20:
574* IF (wantz .and. INDEIG ) Z( 1,1) = Z(1,1) / 1.002 + .002
575*
576*
577 work( 1 ) = lwmin
578 iwork( 1 ) = liwmin
579 RETURN
580*
581* End of DSTEVR
582*
583 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
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 dstevr(jobz, range, n, d, e, vl, vu, il, iu, abstol, m, w, z, ldz, isuppz, work, lwork, iwork, liwork, info)
DSTEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition dstevr.f:303
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82