LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
zhbgvx.f
Go to the documentation of this file.
1 *> \brief \b ZHBGVX
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download ZHBGVX + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zhbgvx.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zhbgvx.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zhbgvx.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE ZHBGVX( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
22 * LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z,
23 * LDZ, WORK, RWORK, IWORK, IFAIL, INFO )
24 *
25 * .. Scalar Arguments ..
26 * CHARACTER JOBZ, RANGE, UPLO
27 * INTEGER IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
28 * $ N
29 * DOUBLE PRECISION ABSTOL, VL, VU
30 * ..
31 * .. Array Arguments ..
32 * INTEGER IFAIL( * ), IWORK( * )
33 * DOUBLE PRECISION RWORK( * ), W( * )
34 * COMPLEX*16 AB( LDAB, * ), BB( LDBB, * ), Q( LDQ, * ),
35 * $ WORK( * ), Z( LDZ, * )
36 * ..
37 *
38 *
39 *> \par Purpose:
40 * =============
41 *>
42 *> \verbatim
43 *>
44 *> ZHBGVX computes all the eigenvalues, and optionally, the eigenvectors
45 *> of a complex generalized Hermitian-definite banded eigenproblem, of
46 *> the form A*x=(lambda)*B*x. Here A and B are assumed to be Hermitian
47 *> and banded, and B is also positive definite. Eigenvalues and
48 *> eigenvectors can be selected by specifying either all eigenvalues,
49 *> a range of values or a range of indices for the desired eigenvalues.
50 *> \endverbatim
51 *
52 * Arguments:
53 * ==========
54 *
55 *> \param[in] JOBZ
56 *> \verbatim
57 *> JOBZ is CHARACTER*1
58 *> = 'N': Compute eigenvalues only;
59 *> = 'V': Compute eigenvalues and eigenvectors.
60 *> \endverbatim
61 *>
62 *> \param[in] RANGE
63 *> \verbatim
64 *> RANGE is CHARACTER*1
65 *> = 'A': all eigenvalues will be found;
66 *> = 'V': all eigenvalues in the half-open interval (VL,VU]
67 *> will be found;
68 *> = 'I': the IL-th through IU-th eigenvalues will be found.
69 *> \endverbatim
70 *>
71 *> \param[in] UPLO
72 *> \verbatim
73 *> UPLO is CHARACTER*1
74 *> = 'U': Upper triangles of A and B are stored;
75 *> = 'L': Lower triangles of A and B are stored.
76 *> \endverbatim
77 *>
78 *> \param[in] N
79 *> \verbatim
80 *> N is INTEGER
81 *> The order of the matrices A and B. N >= 0.
82 *> \endverbatim
83 *>
84 *> \param[in] KA
85 *> \verbatim
86 *> KA is INTEGER
87 *> The number of superdiagonals of the matrix A if UPLO = 'U',
88 *> or the number of subdiagonals if UPLO = 'L'. KA >= 0.
89 *> \endverbatim
90 *>
91 *> \param[in] KB
92 *> \verbatim
93 *> KB is INTEGER
94 *> The number of superdiagonals of the matrix B if UPLO = 'U',
95 *> or the number of subdiagonals if UPLO = 'L'. KB >= 0.
96 *> \endverbatim
97 *>
98 *> \param[in,out] AB
99 *> \verbatim
100 *> AB is COMPLEX*16 array, dimension (LDAB, N)
101 *> On entry, the upper or lower triangle of the Hermitian band
102 *> matrix A, stored in the first ka+1 rows of the array. The
103 *> j-th column of A is stored in the j-th column of the array AB
104 *> as follows:
105 *> if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
106 *> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+ka).
107 *>
108 *> On exit, the contents of AB are destroyed.
109 *> \endverbatim
110 *>
111 *> \param[in] LDAB
112 *> \verbatim
113 *> LDAB is INTEGER
114 *> The leading dimension of the array AB. LDAB >= KA+1.
115 *> \endverbatim
116 *>
117 *> \param[in,out] BB
118 *> \verbatim
119 *> BB is COMPLEX*16 array, dimension (LDBB, N)
120 *> On entry, the upper or lower triangle of the Hermitian band
121 *> matrix B, stored in the first kb+1 rows of the array. The
122 *> j-th column of B is stored in the j-th column of the array BB
123 *> as follows:
124 *> if UPLO = 'U', BB(kb+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
125 *> if UPLO = 'L', BB(1+i-j,j) = B(i,j) for j<=i<=min(n,j+kb).
126 *>
127 *> On exit, the factor S from the split Cholesky factorization
128 *> B = S**H*S, as returned by ZPBSTF.
129 *> \endverbatim
130 *>
131 *> \param[in] LDBB
132 *> \verbatim
133 *> LDBB is INTEGER
134 *> The leading dimension of the array BB. LDBB >= KB+1.
135 *> \endverbatim
136 *>
137 *> \param[out] Q
138 *> \verbatim
139 *> Q is COMPLEX*16 array, dimension (LDQ, N)
140 *> If JOBZ = 'V', the n-by-n matrix used in the reduction of
141 *> A*x = (lambda)*B*x to standard form, i.e. C*x = (lambda)*x,
142 *> and consequently C to tridiagonal form.
143 *> If JOBZ = 'N', the array Q is not referenced.
144 *> \endverbatim
145 *>
146 *> \param[in] LDQ
147 *> \verbatim
148 *> LDQ is INTEGER
149 *> The leading dimension of the array Q. If JOBZ = 'N',
150 *> LDQ >= 1. If JOBZ = 'V', LDQ >= max(1,N).
151 *> \endverbatim
152 *>
153 *> \param[in] VL
154 *> \verbatim
155 *> VL is DOUBLE PRECISION
156 *>
157 *> If RANGE='V', the lower bound of the interval to
158 *> be searched for eigenvalues. VL < VU.
159 *> Not referenced if RANGE = 'A' or 'I'.
160 *> \endverbatim
161 *>
162 *> \param[in] VU
163 *> \verbatim
164 *> VU is DOUBLE PRECISION
165 *>
166 *> If RANGE='V', the upper bound of the interval to
167 *> be searched for eigenvalues. VL < VU.
168 *> Not referenced if RANGE = 'A' or 'I'.
169 *> \endverbatim
170 *>
171 *> \param[in] IL
172 *> \verbatim
173 *> IL is INTEGER
174 *>
175 *> If RANGE='I', the index of the
176 *> smallest eigenvalue to be returned.
177 *> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
178 *> Not referenced if RANGE = 'A' or 'V'.
179 *> \endverbatim
180 *>
181 *> \param[in] IU
182 *> \verbatim
183 *> IU is INTEGER
184 *>
185 *> If RANGE='I', the index of the
186 *> largest eigenvalue to be returned.
187 *> 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
188 *> Not referenced if RANGE = 'A' or 'V'.
189 *> \endverbatim
190 *>
191 *> \param[in] ABSTOL
192 *> \verbatim
193 *> ABSTOL is DOUBLE PRECISION
194 *> The absolute error tolerance for the eigenvalues.
195 *> An approximate eigenvalue is accepted as converged
196 *> when it is determined to lie in an interval [a,b]
197 *> of width less than or equal to
198 *>
199 *> ABSTOL + EPS * max( |a|,|b| ) ,
200 *>
201 *> where EPS is the machine precision. If ABSTOL is less than
202 *> or equal to zero, then EPS*|T| will be used in its place,
203 *> where |T| is the 1-norm of the tridiagonal matrix obtained
204 *> by reducing AP to tridiagonal form.
205 *>
206 *> Eigenvalues will be computed most accurately when ABSTOL is
207 *> set to twice the underflow threshold 2*DLAMCH('S'), not zero.
208 *> If this routine returns with INFO>0, indicating that some
209 *> eigenvectors did not converge, try setting ABSTOL to
210 *> 2*DLAMCH('S').
211 *> \endverbatim
212 *>
213 *> \param[out] M
214 *> \verbatim
215 *> M is INTEGER
216 *> The total number of eigenvalues found. 0 <= M <= N.
217 *> If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
218 *> \endverbatim
219 *>
220 *> \param[out] W
221 *> \verbatim
222 *> W is DOUBLE PRECISION array, dimension (N)
223 *> If INFO = 0, the eigenvalues in ascending order.
224 *> \endverbatim
225 *>
226 *> \param[out] Z
227 *> \verbatim
228 *> Z is COMPLEX*16 array, dimension (LDZ, N)
229 *> If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
230 *> eigenvectors, with the i-th column of Z holding the
231 *> eigenvector associated with W(i). The eigenvectors are
232 *> normalized so that Z**H*B*Z = I.
233 *> If JOBZ = 'N', then Z is not referenced.
234 *> \endverbatim
235 *>
236 *> \param[in] LDZ
237 *> \verbatim
238 *> LDZ is INTEGER
239 *> The leading dimension of the array Z. LDZ >= 1, and if
240 *> JOBZ = 'V', LDZ >= N.
241 *> \endverbatim
242 *>
243 *> \param[out] WORK
244 *> \verbatim
245 *> WORK is COMPLEX*16 array, dimension (N)
246 *> \endverbatim
247 *>
248 *> \param[out] RWORK
249 *> \verbatim
250 *> RWORK is DOUBLE PRECISION array, dimension (7*N)
251 *> \endverbatim
252 *>
253 *> \param[out] IWORK
254 *> \verbatim
255 *> IWORK is INTEGER array, dimension (5*N)
256 *> \endverbatim
257 *>
258 *> \param[out] IFAIL
259 *> \verbatim
260 *> IFAIL is INTEGER array, dimension (N)
261 *> If JOBZ = 'V', then if INFO = 0, the first M elements of
262 *> IFAIL are zero. If INFO > 0, then IFAIL contains the
263 *> indices of the eigenvectors that failed to converge.
264 *> If JOBZ = 'N', then IFAIL is not referenced.
265 *> \endverbatim
266 *>
267 *> \param[out] INFO
268 *> \verbatim
269 *> INFO is INTEGER
270 *> = 0: successful exit
271 *> < 0: if INFO = -i, the i-th argument had an illegal value
272 *> > 0: if INFO = i, and i is:
273 *> <= N: then i eigenvectors failed to converge. Their
274 *> indices are stored in array IFAIL.
275 *> > N: if INFO = N + i, for 1 <= i <= N, then ZPBSTF
276 *> returned INFO = i: B is not positive definite.
277 *> The factorization of B could not be completed and
278 *> no eigenvalues or eigenvectors were computed.
279 *> \endverbatim
280 *
281 * Authors:
282 * ========
283 *
284 *> \author Univ. of Tennessee
285 *> \author Univ. of California Berkeley
286 *> \author Univ. of Colorado Denver
287 *> \author NAG Ltd.
288 *
289 *> \date June 2016
290 *
291 *> \ingroup complex16OTHEReigen
292 *
293 *> \par Contributors:
294 * ==================
295 *>
296 *> Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
297 *
298 * =====================================================================
299  SUBROUTINE zhbgvx( JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB,
300  $ ldbb, q, ldq, vl, vu, il, iu, abstol, m, w, z,
301  $ ldz, work, rwork, iwork, ifail, info )
302 *
303 * -- LAPACK driver routine (version 3.6.1) --
304 * -- LAPACK is a software package provided by Univ. of Tennessee, --
305 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
306 * June 2016
307 *
308 * .. Scalar Arguments ..
309  CHARACTER JOBZ, RANGE, UPLO
310  INTEGER IL, INFO, IU, KA, KB, LDAB, LDBB, LDQ, LDZ, M,
311  $ n
312  DOUBLE PRECISION ABSTOL, VL, VU
313 * ..
314 * .. Array Arguments ..
315  INTEGER IFAIL( * ), IWORK( * )
316  DOUBLE PRECISION RWORK( * ), W( * )
317  COMPLEX*16 AB( ldab, * ), BB( ldbb, * ), Q( ldq, * ),
318  $ work( * ), z( ldz, * )
319 * ..
320 *
321 * =====================================================================
322 *
323 * .. Parameters ..
324  DOUBLE PRECISION ZERO
325  parameter ( zero = 0.0d+0 )
326  COMPLEX*16 CZERO, CONE
327  parameter ( czero = ( 0.0d+0, 0.0d+0 ),
328  $ cone = ( 1.0d+0, 0.0d+0 ) )
329 * ..
330 * .. Local Scalars ..
331  LOGICAL ALLEIG, INDEIG, TEST, UPPER, VALEIG, WANTZ
332  CHARACTER ORDER, VECT
333  INTEGER I, IINFO, INDD, INDE, INDEE, INDIBL, INDISP,
334  $ indiwk, indrwk, indwrk, itmp1, j, jj, nsplit
335  DOUBLE PRECISION TMP1
336 * ..
337 * .. External Functions ..
338  LOGICAL LSAME
339  EXTERNAL lsame
340 * ..
341 * .. External Subroutines ..
342  EXTERNAL dcopy, dstebz, dsterf, xerbla, zcopy, zgemv,
344  $ zswap
345 * ..
346 * .. Intrinsic Functions ..
347  INTRINSIC min
348 * ..
349 * .. Executable Statements ..
350 *
351 * Test the input parameters.
352 *
353  wantz = lsame( jobz, 'V' )
354  upper = lsame( uplo, 'U' )
355  alleig = lsame( range, 'A' )
356  valeig = lsame( range, 'V' )
357  indeig = lsame( range, 'I' )
358 *
359  info = 0
360  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
361  info = -1
362  ELSE IF( .NOT.( alleig .OR. valeig .OR. indeig ) ) THEN
363  info = -2
364  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
365  info = -3
366  ELSE IF( n.LT.0 ) THEN
367  info = -4
368  ELSE IF( ka.LT.0 ) THEN
369  info = -5
370  ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
371  info = -6
372  ELSE IF( ldab.LT.ka+1 ) THEN
373  info = -8
374  ELSE IF( ldbb.LT.kb+1 ) THEN
375  info = -10
376  ELSE IF( ldq.LT.1 .OR. ( wantz .AND. ldq.LT.n ) ) THEN
377  info = -12
378  ELSE
379  IF( valeig ) THEN
380  IF( n.GT.0 .AND. vu.LE.vl )
381  $ info = -14
382  ELSE IF( indeig ) THEN
383  IF( il.LT.1 .OR. il.GT.max( 1, n ) ) THEN
384  info = -15
385  ELSE IF ( iu.LT.min( n, il ) .OR. iu.GT.n ) THEN
386  info = -16
387  END IF
388  END IF
389  END IF
390  IF( info.EQ.0) THEN
391  IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
392  info = -21
393  END IF
394  END IF
395 *
396  IF( info.NE.0 ) THEN
397  CALL xerbla( 'ZHBGVX', -info )
398  RETURN
399  END IF
400 *
401 * Quick return if possible
402 *
403  m = 0
404  IF( n.EQ.0 )
405  $ RETURN
406 *
407 * Form a split Cholesky factorization of B.
408 *
409  CALL zpbstf( uplo, n, kb, bb, ldbb, info )
410  IF( info.NE.0 ) THEN
411  info = n + info
412  RETURN
413  END IF
414 *
415 * Transform problem to standard eigenvalue problem.
416 *
417  CALL zhbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, q, ldq,
418  $ work, rwork, iinfo )
419 *
420 * Solve the standard eigenvalue problem.
421 * Reduce Hermitian band matrix to tridiagonal form.
422 *
423  indd = 1
424  inde = indd + n
425  indrwk = inde + n
426  indwrk = 1
427  IF( wantz ) THEN
428  vect = 'U'
429  ELSE
430  vect = 'N'
431  END IF
432  CALL zhbtrd( vect, uplo, n, ka, ab, ldab, rwork( indd ),
433  $ rwork( inde ), q, ldq, work( indwrk ), iinfo )
434 *
435 * If all eigenvalues are desired and ABSTOL is less than or equal
436 * to zero, then call DSTERF or ZSTEQR. If this fails for some
437 * eigenvalue, then try DSTEBZ.
438 *
439  test = .false.
440  IF( indeig ) THEN
441  IF( il.EQ.1 .AND. iu.EQ.n ) THEN
442  test = .true.
443  END IF
444  END IF
445  IF( ( alleig .OR. test ) .AND. ( abstol.LE.zero ) ) THEN
446  CALL dcopy( n, rwork( indd ), 1, w, 1 )
447  indee = indrwk + 2*n
448  CALL dcopy( n-1, rwork( inde ), 1, rwork( indee ), 1 )
449  IF( .NOT.wantz ) THEN
450  CALL dsterf( n, w, rwork( indee ), info )
451  ELSE
452  CALL zlacpy( 'A', n, n, q, ldq, z, ldz )
453  CALL zsteqr( jobz, n, w, rwork( indee ), z, ldz,
454  $ rwork( indrwk ), info )
455  IF( info.EQ.0 ) THEN
456  DO 10 i = 1, n
457  ifail( i ) = 0
458  10 CONTINUE
459  END IF
460  END IF
461  IF( info.EQ.0 ) THEN
462  m = n
463  GO TO 30
464  END IF
465  info = 0
466  END IF
467 *
468 * Otherwise, call DSTEBZ and, if eigenvectors are desired,
469 * call ZSTEIN.
470 *
471  IF( wantz ) THEN
472  order = 'B'
473  ELSE
474  order = 'E'
475  END IF
476  indibl = 1
477  indisp = indibl + n
478  indiwk = indisp + n
479  CALL dstebz( range, order, n, vl, vu, il, iu, abstol,
480  $ rwork( indd ), rwork( inde ), m, nsplit, w,
481  $ iwork( indibl ), iwork( indisp ), rwork( indrwk ),
482  $ iwork( indiwk ), info )
483 *
484  IF( wantz ) THEN
485  CALL zstein( n, rwork( indd ), rwork( inde ), m, w,
486  $ iwork( indibl ), iwork( indisp ), z, ldz,
487  $ rwork( indrwk ), iwork( indiwk ), ifail, info )
488 *
489 * Apply unitary matrix used in reduction to tridiagonal
490 * form to eigenvectors returned by ZSTEIN.
491 *
492  DO 20 j = 1, m
493  CALL zcopy( n, z( 1, j ), 1, work( 1 ), 1 )
494  CALL zgemv( 'N', n, n, cone, q, ldq, work, 1, czero,
495  $ z( 1, j ), 1 )
496  20 CONTINUE
497  END IF
498 *
499  30 CONTINUE
500 *
501 * If eigenvalues are not in order, then sort them, along with
502 * eigenvectors.
503 *
504  IF( wantz ) THEN
505  DO 50 j = 1, m - 1
506  i = 0
507  tmp1 = w( j )
508  DO 40 jj = j + 1, m
509  IF( w( jj ).LT.tmp1 ) THEN
510  i = jj
511  tmp1 = w( jj )
512  END IF
513  40 CONTINUE
514 *
515  IF( i.NE.0 ) THEN
516  itmp1 = iwork( indibl+i-1 )
517  w( i ) = w( j )
518  iwork( indibl+i-1 ) = iwork( indibl+j-1 )
519  w( j ) = tmp1
520  iwork( indibl+j-1 ) = itmp1
521  CALL zswap( n, z( 1, i ), 1, z( 1, j ), 1 )
522  IF( info.NE.0 ) THEN
523  itmp1 = ifail( i )
524  ifail( i ) = ifail( j )
525  ifail( j ) = itmp1
526  END IF
527  END IF
528  50 CONTINUE
529  END IF
530 *
531  RETURN
532 *
533 * End of ZHBGVX
534 *
535  END
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine dstebz(RANGE, ORDER, N, VL, VU, IL, IU, ABSTOL, D, E, M, NSPLIT, W, IBLOCK, ISPLIT, WORK, IWORK, INFO)
DSTEBZ
Definition: dstebz.f:275
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine zcopy(N, ZX, INCX, ZY, INCY)
ZCOPY
Definition: zcopy.f:52
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:160
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:52
subroutine zhbgst(VECT, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, X, LDX, WORK, RWORK, INFO)
ZHBGST
Definition: zhbgst.f:167
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
ZSTEQR
Definition: zsteqr.f:134
subroutine zstein(N, D, E, M, W, IBLOCK, ISPLIT, Z, LDZ, WORK, IWORK, IFAIL, INFO)
ZSTEIN
Definition: zstein.f:184
subroutine zpbstf(UPLO, N, KD, AB, LDAB, INFO)
ZPBSTF
Definition: zpbstf.f:155
subroutine zhbgvx(JOBZ, RANGE, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, Q, LDQ, VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, RWORK, IWORK, IFAIL, INFO)
ZHBGVX
Definition: zhbgvx.f:302
subroutine zhbtrd(VECT, UPLO, N, KD, AB, LDAB, D, E, Q, LDQ, WORK, INFO)
ZHBTRD
Definition: zhbtrd.f:165