LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dsyevx.f
Go to the documentation of this file.
1*> \brief <b> DSYEVX 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 DSYEVX + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsyevx.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsyevx.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsyevx.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DSYEVX( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU,
22* ABSTOL, M, W, Z, LDZ, WORK, LWORK, IWORK,
23* IFAIL, INFO )
24*
25* .. Scalar Arguments ..
26* CHARACTER JOBZ, RANGE, UPLO
27* INTEGER IL, INFO, IU, LDA, LDZ, LWORK, M, N
28* DOUBLE PRECISION ABSTOL, VL, VU
29* ..
30* .. Array Arguments ..
31* INTEGER IFAIL( * ), IWORK( * )
32* DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
33* ..
34*
35*
36*> \par Purpose:
37* =============
38*>
39*> \verbatim
40*>
41*> DSYEVX 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 indices
44*> for the desired eigenvalues.
45*> \endverbatim
46*
47* Arguments:
48* ==========
49*
50*> \param[in] JOBZ
51*> \verbatim
52*> JOBZ is CHARACTER*1
53*> = 'N': Compute eigenvalues only;
54*> = 'V': Compute eigenvalues and eigenvectors.
55*> \endverbatim
56*>
57*> \param[in] RANGE
58*> \verbatim
59*> RANGE is CHARACTER*1
60*> = 'A': all eigenvalues will be found.
61*> = 'V': all eigenvalues in the half-open interval (VL,VU]
62*> will be found.
63*> = 'I': the IL-th through IU-th eigenvalues will be found.
64*> \endverbatim
65*>
66*> \param[in] UPLO
67*> \verbatim
68*> UPLO is CHARACTER*1
69*> = 'U': Upper triangle of A is stored;
70*> = 'L': Lower triangle of A is stored.
71*> \endverbatim
72*>
73*> \param[in] N
74*> \verbatim
75*> N is INTEGER
76*> The order of the matrix A. N >= 0.
77*> \endverbatim
78*>
79*> \param[in,out] A
80*> \verbatim
81*> A is DOUBLE PRECISION array, dimension (LDA, N)
82*> On entry, the symmetric matrix A. If UPLO = 'U', the
83*> leading N-by-N upper triangular part of A contains the
84*> upper triangular part of the matrix A. If UPLO = 'L',
85*> the leading N-by-N lower triangular part of A contains
86*> the lower triangular part of the matrix A.
87*> On exit, the lower triangle (if UPLO='L') or the upper
88*> triangle (if UPLO='U') of A, including the diagonal, is
89*> destroyed.
90*> \endverbatim
91*>
92*> \param[in] LDA
93*> \verbatim
94*> LDA is INTEGER
95*> The leading dimension of the array A. LDA >= max(1,N).
96*> \endverbatim
97*>
98*> \param[in] VL
99*> \verbatim
100*> VL is DOUBLE PRECISION
101*> If RANGE='V', the lower bound of the interval to
102*> be searched for eigenvalues. VL < VU.
103*> Not referenced if RANGE = 'A' or 'I'.
104*> \endverbatim
105*>
106*> \param[in] VU
107*> \verbatim
108*> VU is DOUBLE PRECISION
109*> If RANGE='V', the upper bound of the interval to
110*> be searched for eigenvalues. VL < VU.
111*> Not referenced if RANGE = 'A' or 'I'.
112*> \endverbatim
113*>
114*> \param[in] IL
115*> \verbatim
116*> IL is INTEGER
117*> If RANGE='I', the index of the
118*> smallest eigenvalue to be returned.
119*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
120*> Not referenced if RANGE = 'A' or 'V'.
121*> \endverbatim
122*>
123*> \param[in] IU
124*> \verbatim
125*> IU is INTEGER
126*> If RANGE='I', the index of the
127*> largest eigenvalue to be returned.
128*> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
129*> Not referenced if RANGE = 'A' or 'V'.
130*> \endverbatim
131*>
132*> \param[in] ABSTOL
133*> \verbatim
134*> ABSTOL is DOUBLE PRECISION
135*> The absolute error tolerance for the eigenvalues.
136*> An approximate eigenvalue is accepted as converged
137*> when it is determined to lie in an interval [a,b]
138*> of width less than or equal to
139*>
140*> ABSTOL + EPS * max( |a|,|b| ) ,
141*>
142*> where EPS is the machine precision. If ABSTOL is less than
143*> or equal to zero, then EPS*|T| will be used in its place,
144*> where |T| is the 1-norm of the tridiagonal matrix obtained
145*> by reducing A to tridiagonal form.
146*>
147*> Eigenvalues will be computed most accurately when ABSTOL is
148*> set to twice the underflow threshold 2*DLAMCH('S'), not zero.
149*> If this routine returns with INFO>0, indicating that some
150*> eigenvectors did not converge, try setting ABSTOL to
151*> 2*DLAMCH('S').
152*>
153*> See "Computing Small Singular Values of Bidiagonal Matrices
154*> with Guaranteed High Relative Accuracy," by Demmel and
155*> Kahan, LAPACK Working Note #3.
156*> \endverbatim
157*>
158*> \param[out] M
159*> \verbatim
160*> M is INTEGER
161*> The total number of eigenvalues found. 0 <= M <= N.
162*> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
163*> \endverbatim
164*>
165*> \param[out] W
166*> \verbatim
167*> W is DOUBLE PRECISION array, dimension (N)
168*> On normal exit, the first M elements contain the selected
169*> eigenvalues in ascending order.
170*> \endverbatim
171*>
172*> \param[out] Z
173*> \verbatim
174*> Z is DOUBLE PRECISION array, dimension (LDZ, max(1,M))
175*> If JOBZ = 'V', then if INFO = 0, the first M columns of Z
176*> contain the orthonormal eigenvectors of the matrix A
177*> corresponding to the selected eigenvalues, with the i-th
178*> column of Z holding the eigenvector associated with W(i).
179*> If an eigenvector fails to converge, then that column of Z
180*> contains the latest approximation to the eigenvector, and the
181*> index of the eigenvector is returned in IFAIL.
182*> If JOBZ = 'N', then Z is not referenced.
183*> Note: the user must ensure that at least max(1,M) columns are
184*> supplied in the array Z; if RANGE = 'V', the exact value of M
185*> is not known in advance and an upper bound must be used.
186*> \endverbatim
187*>
188*> \param[in] LDZ
189*> \verbatim
190*> LDZ is INTEGER
191*> The leading dimension of the array Z. LDZ >= 1, and if
192*> JOBZ = 'V', LDZ >= max(1,N).
193*> \endverbatim
194*>
195*> \param[out] WORK
196*> \verbatim
197*> WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
198*> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
199*> \endverbatim
200*>
201*> \param[in] LWORK
202*> \verbatim
203*> LWORK is INTEGER
204*> The length of the array WORK. LWORK >= 1, when N <= 1;
205*> otherwise 8*N.
206*> For optimal efficiency, LWORK >= (NB+3)*N,
207*> where NB is the max of the blocksize for DSYTRD and DORMTR
208*> returned by ILAENV.
209*>
210*> If LWORK = -1, then a workspace query is assumed; the routine
211*> only calculates the optimal size of the WORK array, returns
212*> this value as the first entry of the WORK array, and no error
213*> message related to LWORK is issued by XERBLA.
214*> \endverbatim
215*>
216*> \param[out] IWORK
217*> \verbatim
218*> IWORK is INTEGER array, dimension (5*N)
219*> \endverbatim
220*>
221*> \param[out] IFAIL
222*> \verbatim
223*> IFAIL is INTEGER array, dimension (N)
224*> If JOBZ = 'V', then if INFO = 0, the first M elements of
225*> IFAIL are zero. If INFO > 0, then IFAIL contains the
226*> indices of the eigenvectors that failed to converge.
227*> If JOBZ = 'N', then IFAIL is not referenced.
228*> \endverbatim
229*>
230*> \param[out] INFO
231*> \verbatim
232*> INFO is INTEGER
233*> = 0: successful exit
234*> < 0: if INFO = -i, the i-th argument had an illegal value
235*> > 0: if INFO = i, then i eigenvectors failed to converge.
236*> Their indices are stored in array IFAIL.
237*> \endverbatim
238*
239* Authors:
240* ========
241*
242*> \author Univ. of Tennessee
243*> \author Univ. of California Berkeley
244*> \author Univ. of Colorado Denver
245*> \author NAG Ltd.
246*
247*> \ingroup heevx
248*
249* =====================================================================
250 SUBROUTINE dsyevx( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU,
251 $ ABSTOL, M, W, Z, LDZ, WORK, LWORK, IWORK,
252 $ IFAIL, INFO )
253*
254* -- LAPACK driver routine --
255* -- LAPACK is a software package provided by Univ. of Tennessee, --
256* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
257*
258* .. Scalar Arguments ..
259 CHARACTER JOBZ, RANGE, UPLO
260 INTEGER IL, INFO, IU, LDA, LDZ, LWORK, M, N
261 DOUBLE PRECISION ABSTOL, VL, VU
262* ..
263* .. Array Arguments ..
264 INTEGER IFAIL( * ), IWORK( * )
265 DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ), Z( LDZ, * )
266* ..
267*
268* =====================================================================
269*
270* .. Parameters ..
271 DOUBLE PRECISION ZERO, ONE
272 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
273* ..
274* .. Local Scalars ..
275 LOGICAL ALLEIG, INDEIG, LOWER, LQUERY, TEST, VALEIG,
276 $ WANTZ
277 CHARACTER ORDER
278 INTEGER I, IINFO, IMAX, INDD, INDE, INDEE, INDIBL,
279 $ indisp, indiwo, indtau, indwkn, indwrk, iscale,
280 $ itmp1, j, jj, llwork, llwrkn, lwkmin,
281 $ lwkopt, nb, nsplit
282 DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
283 $ SIGMA, SMLNUM, TMP1, VLL, VUU
284* ..
285* .. External Functions ..
286 LOGICAL LSAME
287 INTEGER ILAENV
288 DOUBLE PRECISION DLAMCH, DLANSY
289 EXTERNAL lsame, ilaenv, dlamch, dlansy
290* ..
291* .. External Subroutines ..
292 EXTERNAL dcopy, dlacpy, dorgtr, dormtr, dscal, dstebz,
294* ..
295* .. Intrinsic Functions ..
296 INTRINSIC max, min, sqrt
297* ..
298* .. Executable Statements ..
299*
300* Test the input parameters.
301*
302 lower = lsame( uplo, 'L' )
303 wantz = lsame( jobz, 'V' )
304 alleig = lsame( range, 'A' )
305 valeig = lsame( range, 'V' )
306 indeig = lsame( range, 'I' )
307 lquery = ( lwork.EQ.-1 )
308*
309 info = 0
310 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
311 info = -1
312 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
313 info = -2
314 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
315 info = -3
316 ELSE IF( n.LT.0 ) THEN
317 info = -4
318 ELSE IF( lda.LT.max( 1, n ) ) THEN
319 info = -6
320 ELSE
321 IF( valeig ) THEN
322 IF( n.GT.0 .AND. vu.LE.vl )
323 $ info = -8
324 ELSE IF( indeig ) THEN
325 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
326 info = -9
327 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
328 info = -10
329 END IF
330 END IF
331 END IF
332 IF( info.EQ.0 ) THEN
333 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
334 info = -15
335 END IF
336 END IF
337*
338 IF( info.EQ.0 ) THEN
339 IF( n.LE.1 ) THEN
340 lwkmin = 1
341 work( 1 ) = lwkmin
342 ELSE
343 lwkmin = 8*n
344 nb = ilaenv( 1, 'DSYTRD', uplo, n, -1, -1, -1 )
345 nb = max( nb, ilaenv( 1, 'DORMTR', uplo, n, -1, -1, -1 ) )
346 lwkopt = max( lwkmin, ( nb + 3 )*n )
347 work( 1 ) = lwkopt
348 END IF
349*
350 IF( lwork.LT.lwkmin .AND. .NOT.lquery )
351 $ info = -17
352 END IF
353*
354 IF( info.NE.0 ) THEN
355 CALL xerbla( 'DSYEVX', -info )
356 RETURN
357 ELSE IF( lquery ) THEN
358 RETURN
359 END IF
360*
361* Quick return if possible
362*
363 m = 0
364 IF( n.EQ.0 ) THEN
365 RETURN
366 END IF
367*
368 IF( n.EQ.1 ) THEN
369 IF( alleig .OR. indeig ) THEN
370 m = 1
371 w( 1 ) = a( 1, 1 )
372 ELSE
373 IF( vl.LT.a( 1, 1 ) .AND. vu.GE.a( 1, 1 ) ) THEN
374 m = 1
375 w( 1 ) = a( 1, 1 )
376 END IF
377 END IF
378 IF( wantz )
379 $ z( 1, 1 ) = one
380 RETURN
381 END IF
382*
383* Get machine constants.
384*
385 safmin = dlamch( 'Safe minimum' )
386 eps = dlamch( 'Precision' )
387 smlnum = safmin / eps
388 bignum = one / smlnum
389 rmin = sqrt( smlnum )
390 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
391*
392* Scale matrix to allowable range, if necessary.
393*
394 iscale = 0
395 abstll = abstol
396 IF( valeig ) THEN
397 vll = vl
398 vuu = vu
399 END IF
400 anrm = dlansy( 'M', uplo, n, a, lda, work )
401 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
402 iscale = 1
403 sigma = rmin / anrm
404 ELSE IF( anrm.GT.rmax ) THEN
405 iscale = 1
406 sigma = rmax / anrm
407 END IF
408 IF( iscale.EQ.1 ) THEN
409 IF( lower ) THEN
410 DO 10 j = 1, n
411 CALL dscal( n-j+1, sigma, a( j, j ), 1 )
412 10 CONTINUE
413 ELSE
414 DO 20 j = 1, n
415 CALL dscal( j, sigma, a( 1, j ), 1 )
416 20 CONTINUE
417 END IF
418 IF( abstol.GT.0 )
419 $ abstll = abstol*sigma
420 IF( valeig ) THEN
421 vll = vl*sigma
422 vuu = vu*sigma
423 END IF
424 END IF
425*
426* Call DSYTRD to reduce symmetric matrix to tridiagonal form.
427*
428 indtau = 1
429 inde = indtau + n
430 indd = inde + n
431 indwrk = indd + n
432 llwork = lwork - indwrk + 1
433 CALL dsytrd( uplo, n, a, lda, work( indd ), work( inde ),
434 $ work( indtau ), work( indwrk ), llwork, iinfo )
435*
436* If all eigenvalues are desired and ABSTOL is less than or equal to
437* zero, then call DSTERF or DORGTR and SSTEQR. If this fails for
438* some eigenvalue, then try DSTEBZ.
439*
440 test = .false.
441 IF( indeig ) THEN
442 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
443 test = .true.
444 END IF
445 END IF
446 IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
447 CALL dcopy( n, work( indd ), 1, w, 1 )
448 indee = indwrk + 2*n
449 IF( .NOT.wantz ) THEN
450 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
451 CALL dsterf( n, w, work( indee ), info )
452 ELSE
453 CALL dlacpy( 'A', n, n, a, lda, z, ldz )
454 CALL dorgtr( uplo, n, z, ldz, work( indtau ),
455 $ work( indwrk ), llwork, iinfo )
456 CALL dcopy( n-1, work( inde ), 1, work( indee ), 1 )
457 CALL dsteqr( jobz, n, w, work( indee ), z, ldz,
458 $ work( indwrk ), info )
459 IF( info.EQ.0 ) THEN
460 DO 30 i = 1, n
461 ifail( i ) = 0
462 30 CONTINUE
463 END IF
464 END IF
465 IF( info.EQ.0 ) THEN
466 m = n
467 GO TO 40
468 END IF
469 info = 0
470 END IF
471*
472* Otherwise, call DSTEBZ and, if eigenvectors are desired, SSTEIN.
473*
474 IF( wantz ) THEN
475 order = 'B'
476 ELSE
477 order = 'E'
478 END IF
479 indibl = 1
480 indisp = indibl + n
481 indiwo = indisp + n
482 CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
483 $ work( indd ), work( inde ), m, nsplit, w,
484 $ iwork( indibl ), iwork( indisp ), work( indwrk ),
485 $ iwork( indiwo ), info )
486*
487 IF( wantz ) THEN
488 CALL dstein( n, work( indd ), work( inde ), m, w,
489 $ iwork( indibl ), iwork( indisp ), z, ldz,
490 $ work( indwrk ), iwork( indiwo ), ifail, info )
491*
492* Apply orthogonal matrix used in reduction to tridiagonal
493* form to eigenvectors returned by DSTEIN.
494*
495 indwkn = inde
496 llwrkn = lwork - indwkn + 1
497 CALL dormtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ), z,
498 $ ldz, work( indwkn ), llwrkn, iinfo )
499 END IF
500*
501* If matrix was scaled, then rescale eigenvalues appropriately.
502*
503 40 CONTINUE
504 IF( iscale.EQ.1 ) THEN
505 IF( info.EQ.0 ) THEN
506 imax = m
507 ELSE
508 imax = info - 1
509 END IF
510 CALL dscal( imax, one / sigma, w, 1 )
511 END IF
512*
513* If eigenvalues are not in order, then sort them, along with
514* eigenvectors.
515*
516 IF( wantz ) THEN
517 DO 60 j = 1, m - 1
518 i = 0
519 tmp1 = w( j )
520 DO 50 jj = j + 1, m
521 IF( w( jj ).LT.tmp1 ) THEN
522 i = jj
523 tmp1 = w( jj )
524 END IF
525 50 CONTINUE
526*
527 IF( i.NE.0 ) THEN
528 itmp1 = iwork( indibl+i-1 )
529 w( i ) = w( j )
530 iwork( indibl+i-1 ) = iwork( indibl+j-1 )
531 w( j ) = tmp1
532 iwork( indibl+j-1 ) = itmp1
533 CALL dswap( n, z( 1, i ), 1, z( 1, j ), 1 )
534 IF( info.NE.0 ) THEN
535 itmp1 = ifail( i )
536 ifail( i ) = ifail( j )
537 ifail( j ) = itmp1
538 END IF
539 END IF
540 60 CONTINUE
541 END IF
542*
543* Set WORK(1) to optimal workspace size.
544*
545 work( 1 ) = lwkopt
546*
547 RETURN
548*
549* End of DSYEVX
550*
551 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine dsyevx(jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz, work, lwork, iwork, ifail, info)
DSYEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices
Definition dsyevx.f:253
subroutine dsytrd(uplo, n, a, lda, d, e, tau, work, lwork, info)
DSYTRD
Definition dsytrd.f:192
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
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 dsteqr(compz, n, d, e, z, ldz, work, info)
DSTEQR
Definition dsteqr.f:131
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82
subroutine dorgtr(uplo, n, a, lda, tau, work, lwork, info)
DORGTR
Definition dorgtr.f:123
subroutine dormtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
DORMTR
Definition dormtr.f:171