LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
chbgvd.f
Go to the documentation of this file.
1*> \brief \b CHBGVD
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download CHBGVD + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/chbgvd.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/chbgvd.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/chbgvd.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE CHBGVD( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W,
22* Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK,
23* LIWORK, INFO )
24*
25* .. Scalar Arguments ..
26* CHARACTER JOBZ, UPLO
27* INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, LIWORK, LRWORK,
28* $ LWORK, N
29* ..
30* .. Array Arguments ..
31* INTEGER IWORK( * )
32* REAL RWORK( * ), W( * )
33* COMPLEX AB( LDAB, * ), BB( LDBB, * ), WORK( * ),
34* $ Z( LDZ, * )
35* ..
36*
37*
38*> \par Purpose:
39* =============
40*>
41*> \verbatim
42*>
43*> CHBGVD computes all the eigenvalues, and optionally, the eigenvectors
44*> of a complex generalized Hermitian-definite banded eigenproblem, of
45*> the form A*x=(lambda)*B*x. Here A and B are assumed to be Hermitian
46*> and banded, and B is also positive definite. If eigenvectors are
47*> desired, it uses a divide and conquer algorithm.
48*>
49*> \endverbatim
50*
51* Arguments:
52* ==========
53*
54*> \param[in] JOBZ
55*> \verbatim
56*> JOBZ is CHARACTER*1
57*> = 'N': Compute eigenvalues only;
58*> = 'V': Compute eigenvalues and eigenvectors.
59*> \endverbatim
60*>
61*> \param[in] UPLO
62*> \verbatim
63*> UPLO is CHARACTER*1
64*> = 'U': Upper triangles of A and B are stored;
65*> = 'L': Lower triangles of A and B are stored.
66*> \endverbatim
67*>
68*> \param[in] N
69*> \verbatim
70*> N is INTEGER
71*> The order of the matrices A and B. N >= 0.
72*> \endverbatim
73*>
74*> \param[in] KA
75*> \verbatim
76*> KA is INTEGER
77*> The number of superdiagonals of the matrix A if UPLO = 'U',
78*> or the number of subdiagonals if UPLO = 'L'. KA >= 0.
79*> \endverbatim
80*>
81*> \param[in] KB
82*> \verbatim
83*> KB is INTEGER
84*> The number of superdiagonals of the matrix B if UPLO = 'U',
85*> or the number of subdiagonals if UPLO = 'L'. KB >= 0.
86*> \endverbatim
87*>
88*> \param[in,out] AB
89*> \verbatim
90*> AB is COMPLEX array, dimension (LDAB, N)
91*> On entry, the upper or lower triangle of the Hermitian band
92*> matrix A, stored in the first ka+1 rows of the array. The
93*> j-th column of A is stored in the j-th column of the array AB
94*> as follows:
95*> if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
96*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+ka).
97*>
98*> On exit, the contents of AB are destroyed.
99*> \endverbatim
100*>
101*> \param[in] LDAB
102*> \verbatim
103*> LDAB is INTEGER
104*> The leading dimension of the array AB. LDAB >= KA+1.
105*> \endverbatim
106*>
107*> \param[in,out] BB
108*> \verbatim
109*> BB is COMPLEX array, dimension (LDBB, N)
110*> On entry, the upper or lower triangle of the Hermitian band
111*> matrix B, stored in the first kb+1 rows of the array. The
112*> j-th column of B is stored in the j-th column of the array BB
113*> as follows:
114*> if UPLO = 'U', BB(kb+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
115*> if UPLO = 'L', BB(1+i-j,j) = B(i,j) for j<=i<=min(n,j+kb).
116*>
117*> On exit, the factor S from the split Cholesky factorization
118*> B = S**H*S, as returned by CPBSTF.
119*> \endverbatim
120*>
121*> \param[in] LDBB
122*> \verbatim
123*> LDBB is INTEGER
124*> The leading dimension of the array BB. LDBB >= KB+1.
125*> \endverbatim
126*>
127*> \param[out] W
128*> \verbatim
129*> W is REAL array, dimension (N)
130*> If INFO = 0, the eigenvalues in ascending order.
131*> \endverbatim
132*>
133*> \param[out] Z
134*> \verbatim
135*> Z is COMPLEX array, dimension (LDZ, N)
136*> If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
137*> eigenvectors, with the i-th column of Z holding the
138*> eigenvector associated with W(i). The eigenvectors are
139*> normalized so that Z**H*B*Z = I.
140*> If JOBZ = 'N', then Z is not referenced.
141*> \endverbatim
142*>
143*> \param[in] LDZ
144*> \verbatim
145*> LDZ is INTEGER
146*> The leading dimension of the array Z. LDZ >= 1, and if
147*> JOBZ = 'V', LDZ >= N.
148*> \endverbatim
149*>
150*> \param[out] WORK
151*> \verbatim
152*> WORK is COMPLEX array, dimension (MAX(1,LWORK))
153*> On exit, if INFO=0, WORK(1) returns the optimal LWORK.
154*> \endverbatim
155*>
156*> \param[in] LWORK
157*> \verbatim
158*> LWORK is INTEGER
159*> The dimension of the array WORK.
160*> If N <= 1, LWORK >= 1.
161*> If JOBZ = 'N' and N > 1, LWORK >= N.
162*> If JOBZ = 'V' and N > 1, LWORK >= 2*N**2.
163*>
164*> If LWORK = -1, then a workspace query is assumed; the routine
165*> only calculates the optimal sizes of the WORK, RWORK and
166*> IWORK arrays, returns these values as the first entries of
167*> the WORK, RWORK and IWORK arrays, and no error message
168*> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
169*> \endverbatim
170*>
171*> \param[out] RWORK
172*> \verbatim
173*> RWORK is REAL array, dimension (MAX(1,LRWORK))
174*> On exit, if INFO=0, RWORK(1) returns the optimal LRWORK.
175*> \endverbatim
176*>
177*> \param[in] LRWORK
178*> \verbatim
179*> LRWORK is INTEGER
180*> The dimension of array RWORK.
181*> If N <= 1, LRWORK >= 1.
182*> If JOBZ = 'N' and N > 1, LRWORK >= N.
183*> If JOBZ = 'V' and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
184*>
185*> If LRWORK = -1, then a workspace query is assumed; the
186*> routine only calculates the optimal sizes of the WORK, RWORK
187*> and IWORK arrays, returns these values as the first entries
188*> of the WORK, RWORK and IWORK arrays, and no error message
189*> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
190*> \endverbatim
191*>
192*> \param[out] IWORK
193*> \verbatim
194*> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
195*> On exit, if INFO=0, IWORK(1) returns the optimal LIWORK.
196*> \endverbatim
197*>
198*> \param[in] LIWORK
199*> \verbatim
200*> LIWORK is INTEGER
201*> The dimension of array IWORK.
202*> If JOBZ = 'N' or N <= 1, LIWORK >= 1.
203*> If JOBZ = 'V' and N > 1, LIWORK >= 3 + 5*N.
204*>
205*> If LIWORK = -1, then a workspace query is assumed; the
206*> routine only calculates the optimal sizes of the WORK, RWORK
207*> and IWORK arrays, returns these values as the first entries
208*> of the WORK, RWORK and IWORK arrays, and no error message
209*> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
210*> \endverbatim
211*>
212*> \param[out] INFO
213*> \verbatim
214*> INFO is INTEGER
215*> = 0: successful exit
216*> < 0: if INFO = -i, the i-th argument had an illegal value
217*> > 0: if INFO = i, and i is:
218*> <= N: the algorithm failed to converge:
219*> i off-diagonal elements of an intermediate
220*> tridiagonal form did not converge to zero;
221*> > N: if INFO = N + i, for 1 <= i <= N, then CPBSTF
222*> returned INFO = i: B is not positive definite.
223*> The factorization of B could not be completed and
224*> no eigenvalues or eigenvectors were computed.
225*> \endverbatim
226*
227* Authors:
228* ========
229*
230*> \author Univ. of Tennessee
231*> \author Univ. of California Berkeley
232*> \author Univ. of Colorado Denver
233*> \author NAG Ltd.
234*
235*> \ingroup hbgvd
236*
237*> \par Contributors:
238* ==================
239*>
240*> Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
241*
242* =====================================================================
243 SUBROUTINE chbgvd( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W,
244 $ Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK,
245 $ LIWORK, INFO )
246*
247* -- LAPACK driver routine --
248* -- LAPACK is a software package provided by Univ. of Tennessee, --
249* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
250*
251* .. Scalar Arguments ..
252 CHARACTER JOBZ, UPLO
253 INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, LIWORK, LRWORK,
254 $ lwork, n
255* ..
256* .. Array Arguments ..
257 INTEGER IWORK( * )
258 REAL RWORK( * ), W( * )
259 COMPLEX AB( LDAB, * ), BB( LDBB, * ), WORK( * ),
260 $ z( ldz, * )
261* ..
262*
263* =====================================================================
264*
265* .. Parameters ..
266 COMPLEX CONE, CZERO
267 PARAMETER ( CONE = ( 1.0e+0, 0.0e+0 ),
268 $ czero = ( 0.0e+0, 0.0e+0 ) )
269* ..
270* .. Local Scalars ..
271 LOGICAL LQUERY, UPPER, WANTZ
272 CHARACTER VECT
273 INTEGER IINFO, INDE, INDWK2, INDWRK, LIWMIN, LLRWK,
274 $ llwk2, lrwmin, lwmin
275* ..
276* .. External Functions ..
277 LOGICAL LSAME
278 REAL SROUNDUP_LWORK
279 EXTERNAL lsame, sroundup_lwork
280* ..
281* .. External Subroutines ..
282 EXTERNAL ssterf, xerbla, cgemm, chbgst, chbtrd, clacpy,
283 $ cpbstf, cstedc
284* ..
285* .. Executable Statements ..
286*
287* Test the input parameters.
288*
289 wantz = lsame( jobz, 'V' )
290 upper = lsame( uplo, 'U' )
291 lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
292*
293 info = 0
294 IF( n.LE.1 ) THEN
295 lwmin = 1+n
296 lrwmin = 1+n
297 liwmin = 1
298 ELSE IF( wantz ) THEN
299 lwmin = 2*n**2
300 lrwmin = 1 + 5*n + 2*n**2
301 liwmin = 3 + 5*n
302 ELSE
303 lwmin = n
304 lrwmin = n
305 liwmin = 1
306 END IF
307 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
308 info = -1
309 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
310 info = -2
311 ELSE IF( n.LT.0 ) THEN
312 info = -3
313 ELSE IF( ka.LT.0 ) THEN
314 info = -4
315 ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
316 info = -5
317 ELSE IF( ldab.LT.ka+1 ) THEN
318 info = -7
319 ELSE IF( ldbb.LT.kb+1 ) THEN
320 info = -9
321 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
322 info = -12
323 END IF
324*
325 IF( info.EQ.0 ) THEN
326 work( 1 ) = sroundup_lwork(lwmin)
327 rwork( 1 ) = lrwmin
328 iwork( 1 ) = liwmin
329*
330 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
331 info = -14
332 ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
333 info = -16
334 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
335 info = -18
336 END IF
337 END IF
338*
339 IF( info.NE.0 ) THEN
340 CALL xerbla( 'CHBGVD', -info )
341 RETURN
342 ELSE IF( lquery ) THEN
343 RETURN
344 END IF
345*
346* Quick return if possible
347*
348 IF( n.EQ.0 )
349 $ RETURN
350*
351* Form a split Cholesky factorization of B.
352*
353 CALL cpbstf( uplo, n, kb, bb, ldbb, info )
354 IF( info.NE.0 ) THEN
355 info = n + info
356 RETURN
357 END IF
358*
359* Transform problem to standard eigenvalue problem.
360*
361 inde = 1
362 indwrk = inde + n
363 indwk2 = 1 + n*n
364 llwk2 = lwork - indwk2 + 2
365 llrwk = lrwork - indwrk + 2
366 CALL chbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, z, ldz,
367 $ work, rwork, iinfo )
368*
369* Reduce Hermitian band matrix to tridiagonal form.
370*
371 IF( wantz ) THEN
372 vect = 'U'
373 ELSE
374 vect = 'N'
375 END IF
376 CALL chbtrd( vect, uplo, n, ka, ab, ldab, w, rwork( inde ), z,
377 $ ldz, work, iinfo )
378*
379* For eigenvalues only, call SSTERF. For eigenvectors, call CSTEDC.
380*
381 IF( .NOT.wantz ) THEN
382 CALL ssterf( n, w, rwork( inde ), info )
383 ELSE
384 CALL cstedc( 'I', n, w, rwork( inde ), work, n, work( indwk2 ),
385 $ llwk2, rwork( indwrk ), llrwk, iwork, liwork,
386 $ info )
387 CALL cgemm( 'N', 'N', n, n, n, cone, z, ldz, work, n, czero,
388 $ work( indwk2 ), n )
389 CALL clacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
390 END IF
391*
392 work( 1 ) = sroundup_lwork(lwmin)
393 rwork( 1 ) = lrwmin
394 iwork( 1 ) = liwmin
395 RETURN
396*
397* End of CHBGVD
398*
399 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
CGEMM
Definition cgemm.f:188
subroutine chbgst(vect, uplo, n, ka, kb, ab, ldab, bb, ldbb, x, ldx, work, rwork, info)
CHBGST
Definition chbgst.f:165
subroutine chbgvd(jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, w, z, ldz, work, lwork, rwork, lrwork, iwork, liwork, info)
CHBGVD
Definition chbgvd.f:246
subroutine chbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
CHBTRD
Definition chbtrd.f:163
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
subroutine cpbstf(uplo, n, kd, ab, ldab, info)
CPBSTF
Definition cpbstf.f:153
subroutine cstedc(compz, n, d, e, z, ldz, work, lwork, rwork, lrwork, iwork, liwork, info)
CSTEDC
Definition cstedc.f:206
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:86