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