LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zheevx.f
Go to the documentation of this file.
1*> \brief <b> ZHEEVX 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 ZHEEVX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zheevx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zheevx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zheevx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZHEEVX( 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* DOUBLE PRECISION ABSTOL, VL, VU
27* ..
28* .. Array Arguments ..
29* INTEGER IFAIL( * ), IWORK( * )
30* DOUBLE PRECISION RWORK( * ), W( * )
31* COMPLEX*16 A( LDA, * ), WORK( * ), Z( LDZ, * )
32* ..
33*
34*
35*> \par Purpose:
36* =============
37*>
38*> \verbatim
39*>
40*> ZHEEVX 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*16 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 DOUBLE PRECISION
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 DOUBLE PRECISION
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 DOUBLE PRECISION
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*DLAMCH('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*DLAMCH('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 DOUBLE PRECISION 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*16 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*16 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 ZHETRD and for
207*> ZUNMTR 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 DOUBLE PRECISION 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 zheevx( 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 DOUBLE PRECISION ABSTOL, VL, VU
267* ..
268* .. Array Arguments ..
269 INTEGER IFAIL( * ), IWORK( * )
270 DOUBLE PRECISION RWORK( * ), W( * )
271 COMPLEX*16 A( LDA, * ), WORK( * ), Z( LDZ, * )
272* ..
273*
274* =====================================================================
275*
276* .. Parameters ..
277 DOUBLE PRECISION ZERO, ONE
278 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
279 COMPLEX*16 CONE
280 parameter( cone = ( 1.0d+0, 0.0d+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 DOUBLE PRECISION ABSTLL, ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN,
291 $ SIGMA, SMLNUM, TMP1, VLL, VUU
292* ..
293* .. External Functions ..
294 LOGICAL LSAME
295 INTEGER ILAENV
296 DOUBLE PRECISION DLAMCH, ZLANHE
297 EXTERNAL lsame, ilaenv, dlamch, zlanhe
298* ..
299* .. External Subroutines ..
300 EXTERNAL dcopy, dscal, dstebz, dsterf, xerbla,
301 $ zdscal,
303 $ zunmtr
304* ..
305* .. Intrinsic Functions ..
306 INTRINSIC dble, max, min, sqrt
307* ..
308* .. Executable Statements ..
309*
310* Test the input parameters.
311*
312 lower = lsame( uplo, 'L' )
313 wantz = lsame( jobz, 'V' )
314 alleig = lsame( range, 'A' )
315 valeig = lsame( range, 'V' )
316 indeig = lsame( range, 'I' )
317 lquery = ( lwork.EQ.-1 )
318*
319 info = 0
320 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
321 info = -1
322 ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
323 info = -2
324 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
325 info = -3
326 ELSE IF( n.LT.0 ) THEN
327 info = -4
328 ELSE IF( lda.LT.max( 1, n ) ) THEN
329 info = -6
330 ELSE
331 IF( valeig ) THEN
332 IF( n.GT.0 .AND. vu.LE.vl )
333 $ info = -8
334 ELSE IF( indeig ) THEN
335 IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
336 info = -9
337 ELSE IF( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
338 info = -10
339 END IF
340 END IF
341 END IF
342 IF( info.EQ.0 ) THEN
343 IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
344 info = -15
345 END IF
346 END IF
347*
348 IF( info.EQ.0 ) THEN
349 IF( n.LE.1 ) THEN
350 lwkmin = 1
351 work( 1 ) = lwkmin
352 ELSE
353 lwkmin = 2*n
354 nb = ilaenv( 1, 'ZHETRD', uplo, n, -1, -1, -1 )
355 nb = max( nb, ilaenv( 1, 'ZUNMTR', uplo, n, -1, -1,
356 $ -1 ) )
357 lwkopt = max( 1, ( nb + 1 )*n )
358 work( 1 ) = lwkopt
359 END IF
360*
361 IF( lwork.LT.lwkmin .AND. .NOT.lquery )
362 $ info = -17
363 END IF
364*
365 IF( info.NE.0 ) THEN
366 CALL xerbla( 'ZHEEVX', -info )
367 RETURN
368 ELSE IF( lquery ) THEN
369 RETURN
370 END IF
371*
372* Quick return if possible
373*
374 m = 0
375 IF( n.EQ.0 ) THEN
376 RETURN
377 END IF
378*
379 IF( n.EQ.1 ) THEN
380 IF( alleig .OR. indeig ) THEN
381 m = 1
382 w( 1 ) = dble( a( 1, 1 ) )
383 ELSE IF( valeig ) THEN
384 IF( vl.LT.dble( a( 1, 1 ) ) .AND. vu.GE.dble( a( 1, 1 ) ) )
385 $ THEN
386 m = 1
387 w( 1 ) = dble( a( 1, 1 ) )
388 END IF
389 END IF
390 IF( wantz )
391 $ z( 1, 1 ) = cone
392 RETURN
393 END IF
394*
395* Get machine constants.
396*
397 safmin = dlamch( 'Safe minimum' )
398 eps = dlamch( 'Precision' )
399 smlnum = safmin / eps
400 bignum = one / smlnum
401 rmin = sqrt( smlnum )
402 rmax = min( sqrt( bignum ), one / sqrt( sqrt( safmin ) ) )
403*
404* Scale matrix to allowable range, if necessary.
405*
406 iscale = 0
407 abstll = abstol
408 IF( valeig ) THEN
409 vll = vl
410 vuu = vu
411 END IF
412 anrm = zlanhe( 'M', uplo, n, a, lda, rwork )
413 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
414 iscale = 1
415 sigma = rmin / anrm
416 ELSE IF( anrm.GT.rmax ) THEN
417 iscale = 1
418 sigma = rmax / anrm
419 END IF
420 IF( iscale.EQ.1 ) THEN
421 IF( lower ) THEN
422 DO 10 j = 1, n
423 CALL zdscal( n-j+1, sigma, a( j, j ), 1 )
424 10 CONTINUE
425 ELSE
426 DO 20 j = 1, n
427 CALL zdscal( j, sigma, a( 1, j ), 1 )
428 20 CONTINUE
429 END IF
430 IF( abstol.GT.0 )
431 $ abstll = abstol*sigma
432 IF( valeig ) THEN
433 vll = vl*sigma
434 vuu = vu*sigma
435 END IF
436 END IF
437*
438* Call ZHETRD to reduce Hermitian matrix to tridiagonal form.
439*
440 indd = 1
441 inde = indd + n
442 indrwk = inde + n
443 indtau = 1
444 indwrk = indtau + n
445 llwork = lwork - indwrk + 1
446 CALL zhetrd( uplo, n, a, lda, rwork( indd ), rwork( inde ),
447 $ work( indtau ), work( indwrk ), llwork, iinfo )
448*
449* If all eigenvalues are desired and ABSTOL is less than or equal to
450* zero, then call DSTERF or ZUNGTR and ZSTEQR. If this fails for
451* some eigenvalue, then try DSTEBZ.
452*
453 test = .false.
454 IF( indeig ) THEN
455 IF( il.EQ.1 .AND. iu.EQ.n ) THEN
456 test = .true.
457 END IF
458 END IF
459 IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
460 CALL dcopy( n, rwork( indd ), 1, w, 1 )
461 indee = indrwk + 2*n
462 IF( .NOT.wantz ) THEN
463 CALL dcopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
464 CALL dsterf( n, w, rwork( indee ), info )
465 ELSE
466 CALL zlacpy( 'A', n, n, a, lda, z, ldz )
467 CALL zungtr( uplo, n, z, ldz, work( indtau ),
468 $ work( indwrk ), llwork, iinfo )
469 CALL dcopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
470 CALL zsteqr( jobz, n, w, rwork( indee ), z, ldz,
471 $ rwork( indrwk ), info )
472 IF( info.EQ.0 ) THEN
473 DO 30 i = 1, n
474 ifail( i ) = 0
475 30 CONTINUE
476 END IF
477 END IF
478 IF( info.EQ.0 ) THEN
479 m = n
480 GO TO 40
481 END IF
482 info = 0
483 END IF
484*
485* Otherwise, call DSTEBZ and, if eigenvectors are desired, ZSTEIN.
486*
487 IF( wantz ) THEN
488 order = 'B'
489 ELSE
490 order = 'E'
491 END IF
492 indibl = 1
493 indisp = indibl + n
494 indiwk = indisp + n
495 CALL dstebz( range, order, n, vll, vuu, il, iu, abstll,
496 $ rwork( indd ), rwork( inde ), m, nsplit, w,
497 $ iwork( indibl ), iwork( indisp ), rwork( indrwk ),
498 $ iwork( indiwk ), info )
499*
500 IF( wantz ) THEN
501 CALL zstein( n, rwork( indd ), rwork( inde ), m, w,
502 $ iwork( indibl ), iwork( indisp ), z, ldz,
503 $ rwork( indrwk ), iwork( indiwk ), ifail, info )
504*
505* Apply unitary matrix used in reduction to tridiagonal
506* form to eigenvectors returned by ZSTEIN.
507*
508 CALL zunmtr( 'L', uplo, 'N', n, m, a, lda, work( indtau ),
509 $ z,
510 $ ldz, work( indwrk ), llwork, iinfo )
511 END IF
512*
513* If matrix was scaled, then rescale eigenvalues appropriately.
514*
515 40 CONTINUE
516 IF( iscale.EQ.1 ) THEN
517 IF( info.EQ.0 ) THEN
518 imax = m
519 ELSE
520 imax = info - 1
521 END IF
522 CALL dscal( imax, one / sigma, w, 1 )
523 END IF
524*
525* If eigenvalues are not in order, then sort them, along with
526* eigenvectors.
527*
528 IF( wantz ) THEN
529 DO 60 j = 1, m - 1
530 i = 0
531 tmp1 = w( j )
532 DO 50 jj = j + 1, m
533 IF( w( jj ).LT.tmp1 ) THEN
534 i = jj
535 tmp1 = w( jj )
536 END IF
537 50 CONTINUE
538*
539 IF( i.NE.0 ) THEN
540 itmp1 = iwork( indibl+i-1 )
541 w( i ) = w( j )
542 iwork( indibl+i-1 ) = iwork( indibl+j-1 )
543 w( j ) = tmp1
544 iwork( indibl+j-1 ) = itmp1
545 CALL zswap( n, z( 1, i ), 1, z( 1, j ), 1 )
546 IF( info.NE.0 ) THEN
547 itmp1 = ifail( i )
548 ifail( i ) = ifail( j )
549 ifail( j ) = itmp1
550 END IF
551 END IF
552 60 CONTINUE
553 END IF
554*
555* Set WORK(1) to optimal complex workspace size.
556*
557 work( 1 ) = lwkopt
558*
559 RETURN
560*
561* End of ZHEEVX
562*
563 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine zheevx(jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz, work, lwork, rwork, iwork, ifail, info)
ZHEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE matrices
Definition zheevx.f:258
subroutine zhetrd(uplo, n, a, lda, d, e, tau, work, lwork, info)
ZHETRD
Definition zhetrd.f:191
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:101
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine zdscal(n, da, zx, incx)
ZDSCAL
Definition zdscal.f:78
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 zstein(n, d, e, m, w, iblock, isplit, z, ldz, work, iwork, ifail, info)
ZSTEIN
Definition zstein.f:180
subroutine zsteqr(compz, n, d, e, z, ldz, work, info)
ZSTEQR
Definition zsteqr.f:130
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:84
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81
subroutine zungtr(uplo, n, a, lda, tau, work, lwork, info)
ZUNGTR
Definition zungtr.f:121
subroutine zunmtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
ZUNMTR
Definition zunmtr.f:170