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