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