LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
chbgv.f
Go to the documentation of this file.
1*> \brief \b CHBGV
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download CHBGV + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/chbgv.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/chbgv.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/chbgv.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE CHBGV( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W, Z,
22* LDZ, WORK, RWORK, INFO )
23*
24* .. Scalar Arguments ..
25* CHARACTER JOBZ, UPLO
26* INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, N
27* ..
28* .. Array Arguments ..
29* REAL RWORK( * ), W( * )
30* COMPLEX AB( LDAB, * ), BB( LDBB, * ), WORK( * ),
31* $ Z( LDZ, * )
32* ..
33*
34*
35*> \par Purpose:
36* =============
37*>
38*> \verbatim
39*>
40*> CHBGV computes all the eigenvalues, and optionally, the eigenvectors
41*> of a complex generalized Hermitian-definite banded eigenproblem, of
42*> the form A*x=(lambda)*B*x. Here A and B are assumed to be Hermitian
43*> and banded, and B is also positive definite.
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] UPLO
57*> \verbatim
58*> UPLO is CHARACTER*1
59*> = 'U': Upper triangles of A and B are stored;
60*> = 'L': Lower triangles of A and B are stored.
61*> \endverbatim
62*>
63*> \param[in] N
64*> \verbatim
65*> N is INTEGER
66*> The order of the matrices A and B. N >= 0.
67*> \endverbatim
68*>
69*> \param[in] KA
70*> \verbatim
71*> KA is INTEGER
72*> The number of superdiagonals of the matrix A if UPLO = 'U',
73*> or the number of subdiagonals if UPLO = 'L'. KA >= 0.
74*> \endverbatim
75*>
76*> \param[in] KB
77*> \verbatim
78*> KB is INTEGER
79*> The number of superdiagonals of the matrix B if UPLO = 'U',
80*> or the number of subdiagonals if UPLO = 'L'. KB >= 0.
81*> \endverbatim
82*>
83*> \param[in,out] AB
84*> \verbatim
85*> AB is COMPLEX array, dimension (LDAB, N)
86*> On entry, the upper or lower triangle of the Hermitian band
87*> matrix A, stored in the first ka+1 rows of the array. The
88*> j-th column of A is stored in the j-th column of the array AB
89*> as follows:
90*> if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
91*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+ka).
92*>
93*> On exit, the contents of AB are destroyed.
94*> \endverbatim
95*>
96*> \param[in] LDAB
97*> \verbatim
98*> LDAB is INTEGER
99*> The leading dimension of the array AB. LDAB >= KA+1.
100*> \endverbatim
101*>
102*> \param[in,out] BB
103*> \verbatim
104*> BB is COMPLEX array, dimension (LDBB, N)
105*> On entry, the upper or lower triangle of the Hermitian band
106*> matrix B, stored in the first kb+1 rows of the array. The
107*> j-th column of B is stored in the j-th column of the array BB
108*> as follows:
109*> if UPLO = 'U', BB(kb+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
110*> if UPLO = 'L', BB(1+i-j,j) = B(i,j) for j<=i<=min(n,j+kb).
111*>
112*> On exit, the factor S from the split Cholesky factorization
113*> B = S**H*S, as returned by CPBSTF.
114*> \endverbatim
115*>
116*> \param[in] LDBB
117*> \verbatim
118*> LDBB is INTEGER
119*> The leading dimension of the array BB. LDBB >= KB+1.
120*> \endverbatim
121*>
122*> \param[out] W
123*> \verbatim
124*> W is REAL array, dimension (N)
125*> If INFO = 0, the eigenvalues in ascending order.
126*> \endverbatim
127*>
128*> \param[out] Z
129*> \verbatim
130*> Z is COMPLEX array, dimension (LDZ, N)
131*> If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
132*> eigenvectors, with the i-th column of Z holding the
133*> eigenvector associated with W(i). The eigenvectors are
134*> normalized so that Z**H*B*Z = I.
135*> If JOBZ = 'N', then Z is not referenced.
136*> \endverbatim
137*>
138*> \param[in] LDZ
139*> \verbatim
140*> LDZ is INTEGER
141*> The leading dimension of the array Z. LDZ >= 1, and if
142*> JOBZ = 'V', LDZ >= N.
143*> \endverbatim
144*>
145*> \param[out] WORK
146*> \verbatim
147*> WORK is COMPLEX array, dimension (N)
148*> \endverbatim
149*>
150*> \param[out] RWORK
151*> \verbatim
152*> RWORK is REAL array, dimension (3*N)
153*> \endverbatim
154*>
155*> \param[out] INFO
156*> \verbatim
157*> INFO is INTEGER
158*> = 0: successful exit
159*> < 0: if INFO = -i, the i-th argument had an illegal value
160*> > 0: if INFO = i, and i is:
161*> <= N: the algorithm failed to converge:
162*> i off-diagonal elements of an intermediate
163*> tridiagonal form did not converge to zero;
164*> > N: if INFO = N + i, for 1 <= i <= N, then CPBSTF
165*> returned INFO = i: B is not positive definite.
166*> The factorization of B could not be completed and
167*> no eigenvalues or eigenvectors were computed.
168*> \endverbatim
169*
170* Authors:
171* ========
172*
173*> \author Univ. of Tennessee
174*> \author Univ. of California Berkeley
175*> \author Univ. of Colorado Denver
176*> \author NAG Ltd.
177*
178*> \ingroup hbgv
179*
180* =====================================================================
181 SUBROUTINE chbgv( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W, Z,
182 $ LDZ, WORK, RWORK, INFO )
183*
184* -- LAPACK driver routine --
185* -- LAPACK is a software package provided by Univ. of Tennessee, --
186* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
187*
188* .. Scalar Arguments ..
189 CHARACTER JOBZ, UPLO
190 INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, N
191* ..
192* .. Array Arguments ..
193 REAL RWORK( * ), W( * )
194 COMPLEX AB( LDAB, * ), BB( LDBB, * ), WORK( * ),
195 $ z( ldz, * )
196* ..
197*
198* =====================================================================
199*
200* .. Local Scalars ..
201 LOGICAL UPPER, WANTZ
202 CHARACTER VECT
203 INTEGER IINFO, INDE, INDWRK
204* ..
205* .. External Functions ..
206 LOGICAL LSAME
207 EXTERNAL lsame
208* ..
209* .. External Subroutines ..
210 EXTERNAL chbgst, chbtrd, cpbstf, csteqr, ssterf, xerbla
211* ..
212* .. Executable Statements ..
213*
214* Test the input parameters.
215*
216 wantz = lsame( jobz, 'V' )
217 upper = lsame( uplo, 'U' )
218*
219 info = 0
220 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
221 info = -1
222 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
223 info = -2
224 ELSE IF( n.LT.0 ) THEN
225 info = -3
226 ELSE IF( ka.LT.0 ) THEN
227 info = -4
228 ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
229 info = -5
230 ELSE IF( ldab.LT.ka+1 ) THEN
231 info = -7
232 ELSE IF( ldbb.LT.kb+1 ) THEN
233 info = -9
234 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
235 info = -12
236 END IF
237 IF( info.NE.0 ) THEN
238 CALL xerbla( 'CHBGV', -info )
239 RETURN
240 END IF
241*
242* Quick return if possible
243*
244 IF( n.EQ.0 )
245 $ RETURN
246*
247* Form a split Cholesky factorization of B.
248*
249 CALL cpbstf( uplo, n, kb, bb, ldbb, info )
250 IF( info.NE.0 ) THEN
251 info = n + info
252 RETURN
253 END IF
254*
255* Transform problem to standard eigenvalue problem.
256*
257 inde = 1
258 indwrk = inde + n
259 CALL chbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, z, ldz,
260 $ work, rwork( indwrk ), iinfo )
261*
262* Reduce to tridiagonal form.
263*
264 IF( wantz ) THEN
265 vect = 'U'
266 ELSE
267 vect = 'N'
268 END IF
269 CALL chbtrd( vect, uplo, n, ka, ab, ldab, w, rwork( inde ), z,
270 $ ldz, work, iinfo )
271*
272* For eigenvalues only, call SSTERF. For eigenvectors, call CSTEQR.
273*
274 IF( .NOT.wantz ) THEN
275 CALL ssterf( n, w, rwork( inde ), info )
276 ELSE
277 CALL csteqr( jobz, n, w, rwork( inde ), z, ldz,
278 $ rwork( indwrk ), info )
279 END IF
280 RETURN
281*
282* End of CHBGV
283*
284 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine chbgst(vect, uplo, n, ka, kb, ab, ldab, bb, ldbb, x, ldx, work, rwork, info)
CHBGST
Definition chbgst.f:165
subroutine chbgv(jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, w, z, ldz, work, rwork, info)
CHBGV
Definition chbgv.f:183
subroutine chbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
CHBTRD
Definition chbtrd.f:163
subroutine cpbstf(uplo, n, kd, ab, ldab, info)
CPBSTF
Definition cpbstf.f:153
subroutine csteqr(compz, n, d, e, z, ldz, work, info)
CSTEQR
Definition csteqr.f:132
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:86