LAPACK 3.12.1
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*> Download CHBGVD + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/chbgvd.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/chbgvd.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/chbgvd.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE CHBGVD( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, W,
20* Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK,
21* LIWORK, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER JOBZ, UPLO
25* INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, LIWORK, LRWORK,
26* $ LWORK, N
27* ..
28* .. Array Arguments ..
29* INTEGER IWORK( * )
30* REAL RWORK( * ), W( * )
31* COMPLEX AB( LDAB, * ), BB( LDBB, * ), WORK( * ),
32* $ Z( LDZ, * )
33* ..
34*
35*
36*> \par Purpose:
37* =============
38*>
39*> \verbatim
40*>
41*> CHBGVD computes all the eigenvalues, and optionally, the eigenvectors
42*> of a complex generalized Hermitian-definite banded eigenproblem, of
43*> the form A*x=(lambda)*B*x. Here A and B are assumed to be Hermitian
44*> and banded, and B is also positive definite. If eigenvectors are
45*> desired, it uses a divide and conquer algorithm.
46*>
47*> \endverbatim
48*
49* Arguments:
50* ==========
51*
52*> \param[in] JOBZ
53*> \verbatim
54*> JOBZ is CHARACTER*1
55*> = 'N': Compute eigenvalues only;
56*> = 'V': Compute eigenvalues and eigenvectors.
57*> \endverbatim
58*>
59*> \param[in] UPLO
60*> \verbatim
61*> UPLO is CHARACTER*1
62*> = 'U': Upper triangles of A and B are stored;
63*> = 'L': Lower triangles of A and B are stored.
64*> \endverbatim
65*>
66*> \param[in] N
67*> \verbatim
68*> N is INTEGER
69*> The order of the matrices A and B. N >= 0.
70*> \endverbatim
71*>
72*> \param[in] KA
73*> \verbatim
74*> KA is INTEGER
75*> The number of superdiagonals of the matrix A if UPLO = 'U',
76*> or the number of subdiagonals if UPLO = 'L'. KA >= 0.
77*> \endverbatim
78*>
79*> \param[in] KB
80*> \verbatim
81*> KB is INTEGER
82*> The number of superdiagonals of the matrix B if UPLO = 'U',
83*> or the number of subdiagonals if UPLO = 'L'. KB >= 0.
84*> \endverbatim
85*>
86*> \param[in,out] AB
87*> \verbatim
88*> AB is COMPLEX array, dimension (LDAB, N)
89*> On entry, the upper or lower triangle of the Hermitian band
90*> matrix A, stored in the first ka+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(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
94*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(n,j+ka).
95*>
96*> On exit, the contents of AB are destroyed.
97*> \endverbatim
98*>
99*> \param[in] LDAB
100*> \verbatim
101*> LDAB is INTEGER
102*> The leading dimension of the array AB. LDAB >= KA+1.
103*> \endverbatim
104*>
105*> \param[in,out] BB
106*> \verbatim
107*> BB is COMPLEX array, dimension (LDBB, N)
108*> On entry, the upper or lower triangle of the Hermitian band
109*> matrix B, stored in the first kb+1 rows of the array. The
110*> j-th column of B is stored in the j-th column of the array BB
111*> as follows:
112*> if UPLO = 'U', BB(kb+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
113*> if UPLO = 'L', BB(1+i-j,j) = B(i,j) for j<=i<=min(n,j+kb).
114*>
115*> On exit, the factor S from the split Cholesky factorization
116*> B = S**H*S, as returned by CPBSTF.
117*> \endverbatim
118*>
119*> \param[in] LDBB
120*> \verbatim
121*> LDBB is INTEGER
122*> The leading dimension of the array BB. LDBB >= KB+1.
123*> \endverbatim
124*>
125*> \param[out] W
126*> \verbatim
127*> W is REAL array, dimension (N)
128*> If INFO = 0, the eigenvalues in ascending order.
129*> \endverbatim
130*>
131*> \param[out] Z
132*> \verbatim
133*> Z is COMPLEX array, dimension (LDZ, N)
134*> If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
135*> eigenvectors, with the i-th column of Z holding the
136*> eigenvector associated with W(i). The eigenvectors are
137*> normalized so that Z**H*B*Z = I.
138*> If JOBZ = 'N', then Z is not referenced.
139*> \endverbatim
140*>
141*> \param[in] LDZ
142*> \verbatim
143*> LDZ is INTEGER
144*> The leading dimension of the array Z. LDZ >= 1, and if
145*> JOBZ = 'V', LDZ >= N.
146*> \endverbatim
147*>
148*> \param[out] WORK
149*> \verbatim
150*> WORK is COMPLEX array, dimension (MAX(1,LWORK))
151*> On exit, if INFO=0, WORK(1) returns the optimal LWORK.
152*> \endverbatim
153*>
154*> \param[in] LWORK
155*> \verbatim
156*> LWORK is INTEGER
157*> The dimension of the array WORK.
158*> If N <= 1, LWORK >= 1.
159*> If JOBZ = 'N' and N > 1, LWORK >= N.
160*> If JOBZ = 'V' and N > 1, LWORK >= 2*N**2.
161*>
162*> If LWORK = -1, then a workspace query is assumed; the routine
163*> only calculates the optimal sizes of the WORK, RWORK and
164*> IWORK arrays, returns these values as the first entries of
165*> the WORK, RWORK and IWORK arrays, and no error message
166*> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
167*> \endverbatim
168*>
169*> \param[out] RWORK
170*> \verbatim
171*> RWORK is REAL array, dimension (MAX(1,LRWORK))
172*> On exit, if INFO=0, RWORK(1) returns the optimal LRWORK.
173*> \endverbatim
174*>
175*> \param[in] LRWORK
176*> \verbatim
177*> LRWORK is INTEGER
178*> The dimension of array RWORK.
179*> If N <= 1, LRWORK >= 1.
180*> If JOBZ = 'N' and N > 1, LRWORK >= N.
181*> If JOBZ = 'V' and N > 1, LRWORK >= 1 + 5*N + 2*N**2.
182*>
183*> If LRWORK = -1, then a workspace query is assumed; the
184*> routine only calculates the optimal sizes of the WORK, RWORK
185*> and IWORK arrays, returns these values as the first entries
186*> of the WORK, RWORK and IWORK arrays, and no error message
187*> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
188*> \endverbatim
189*>
190*> \param[out] IWORK
191*> \verbatim
192*> IWORK is INTEGER array, dimension (MAX(1,LIWORK))
193*> On exit, if INFO=0, IWORK(1) returns the optimal LIWORK.
194*> \endverbatim
195*>
196*> \param[in] LIWORK
197*> \verbatim
198*> LIWORK is INTEGER
199*> The dimension of array IWORK.
200*> If JOBZ = 'N' or N <= 1, LIWORK >= 1.
201*> If JOBZ = 'V' and N > 1, LIWORK >= 3 + 5*N.
202*>
203*> If LIWORK = -1, then a workspace query is assumed; the
204*> routine only calculates the optimal sizes of the WORK, RWORK
205*> and IWORK arrays, returns these values as the first entries
206*> of the WORK, RWORK and IWORK arrays, and no error message
207*> related to LWORK or LRWORK or LIWORK is issued by XERBLA.
208*> \endverbatim
209*>
210*> \param[out] INFO
211*> \verbatim
212*> INFO is INTEGER
213*> = 0: successful exit
214*> < 0: if INFO = -i, the i-th argument had an illegal value
215*> > 0: if INFO = i, and i is:
216*> <= N: the algorithm failed to converge:
217*> i off-diagonal elements of an intermediate
218*> tridiagonal form did not converge to zero;
219*> > N: if INFO = N + i, for 1 <= i <= N, then CPBSTF
220*> returned INFO = i: B is not positive definite.
221*> The factorization of B could not be completed and
222*> no eigenvalues or eigenvectors were computed.
223*> \endverbatim
224*
225* Authors:
226* ========
227*
228*> \author Univ. of Tennessee
229*> \author Univ. of California Berkeley
230*> \author Univ. of Colorado Denver
231*> \author NAG Ltd.
232*
233*> \ingroup hbgvd
234*
235*> \par Contributors:
236* ==================
237*>
238*> Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA
239*
240* =====================================================================
241 SUBROUTINE chbgvd( JOBZ, UPLO, N, KA, KB, AB, LDAB, BB, LDBB,
242 $ W,
243 $ Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK,
244 $ LIWORK, INFO )
245*
246* -- LAPACK driver routine --
247* -- LAPACK is a software package provided by Univ. of Tennessee, --
248* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
249*
250* .. Scalar Arguments ..
251 CHARACTER JOBZ, UPLO
252 INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, LIWORK, LRWORK,
253 $ LWORK, N
254* ..
255* .. Array Arguments ..
256 INTEGER IWORK( * )
257 REAL RWORK( * ), W( * )
258 COMPLEX AB( LDAB, * ), BB( LDBB, * ), WORK( * ),
259 $ z( ldz, * )
260* ..
261*
262* =====================================================================
263*
264* .. Parameters ..
265 COMPLEX CONE, CZERO
266 PARAMETER ( CONE = ( 1.0e+0, 0.0e+0 ),
267 $ czero = ( 0.0e+0, 0.0e+0 ) )
268* ..
269* .. Local Scalars ..
270 LOGICAL LQUERY, UPPER, WANTZ
271 CHARACTER VECT
272 INTEGER IINFO, INDE, INDWK2, INDWRK, LIWMIN, LLRWK,
273 $ llwk2, lrwmin, lwmin
274* ..
275* .. External Functions ..
276 LOGICAL LSAME
277 REAL SROUNDUP_LWORK
278 EXTERNAL LSAME, SROUNDUP_LWORK
279* ..
280* .. External Subroutines ..
281 EXTERNAL ssterf, xerbla, cgemm, chbgst, chbtrd,
282 $ 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 ) = real( 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,
385 $ work( indwk2 ),
386 $ llwk2, rwork( indwrk ), llrwk, iwork, liwork,
387 $ info )
388 CALL cgemm( 'N', 'N', n, n, n, cone, z, ldz, work, n, czero,
389 $ work( indwk2 ), n )
390 CALL clacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
391 END IF
392*
393 work( 1 ) = sroundup_lwork(lwmin)
394 rwork( 1 ) = real( lrwmin )
395 iwork( 1 ) = liwmin
396 RETURN
397*
398* End of CHBGVD
399*
400 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:164
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:245
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 cstedc(compz, n, d, e, z, ldz, work, lwork, rwork, lrwork, iwork, liwork, info)
CSTEDC
Definition cstedc.f:204
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:84