LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ chbgvd()

subroutine chbgvd ( character  jobz,
character  uplo,
integer  n,
integer  ka,
integer  kb,
complex, dimension( ldab, * )  ab,
integer  ldab,
complex, dimension( ldbb, * )  bb,
integer  ldbb,
real, dimension( * )  w,
complex, dimension( ldz, * )  z,
integer  ldz,
complex, dimension( * )  work,
integer  lwork,
real, dimension( * )  rwork,
integer  lrwork,
integer, dimension( * )  iwork,
integer  liwork,
integer  info 
)

CHBGVD

Download CHBGVD + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 CHBGVD computes all the eigenvalues, and optionally, the eigenvectors
 of a complex generalized Hermitian-definite banded eigenproblem, of
 the form A*x=(lambda)*B*x. Here A and B are assumed to be Hermitian
 and banded, and B is also positive definite.  If eigenvectors are
 desired, it uses a divide and conquer algorithm.
Parameters
[in]JOBZ
          JOBZ is CHARACTER*1
          = 'N':  Compute eigenvalues only;
          = 'V':  Compute eigenvalues and eigenvectors.
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangles of A and B are stored;
          = 'L':  Lower triangles of A and B are stored.
[in]N
          N is INTEGER
          The order of the matrices A and B.  N >= 0.
[in]KA
          KA is INTEGER
          The number of superdiagonals of the matrix A if UPLO = 'U',
          or the number of subdiagonals if UPLO = 'L'. KA >= 0.
[in]KB
          KB is INTEGER
          The number of superdiagonals of the matrix B if UPLO = 'U',
          or the number of subdiagonals if UPLO = 'L'. KB >= 0.
[in,out]AB
          AB is COMPLEX array, dimension (LDAB, N)
          On entry, the upper or lower triangle of the Hermitian band
          matrix A, stored in the first ka+1 rows of the array.  The
          j-th column of A is stored in the j-th column of the array AB
          as follows:
          if UPLO = 'U', AB(ka+1+i-j,j) = A(i,j) for max(1,j-ka)<=i<=j;
          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+ka).

          On exit, the contents of AB are destroyed.
[in]LDAB
          LDAB is INTEGER
          The leading dimension of the array AB.  LDAB >= KA+1.
[in,out]BB
          BB is COMPLEX array, dimension (LDBB, N)
          On entry, the upper or lower triangle of the Hermitian band
          matrix B, stored in the first kb+1 rows of the array.  The
          j-th column of B is stored in the j-th column of the array BB
          as follows:
          if UPLO = 'U', BB(kb+1+i-j,j) = B(i,j) for max(1,j-kb)<=i<=j;
          if UPLO = 'L', BB(1+i-j,j)    = B(i,j) for j<=i<=min(n,j+kb).

          On exit, the factor S from the split Cholesky factorization
          B = S**H*S, as returned by CPBSTF.
[in]LDBB
          LDBB is INTEGER
          The leading dimension of the array BB.  LDBB >= KB+1.
[out]W
          W is REAL array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is COMPLEX array, dimension (LDZ, N)
          If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
          eigenvectors, with the i-th column of Z holding the
          eigenvector associated with W(i). The eigenvectors are
          normalized so that Z**H*B*Z = I.
          If JOBZ = 'N', then Z is not referenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= 1, and if
          JOBZ = 'V', LDZ >= N.
[out]WORK
          WORK is COMPLEX array, dimension (MAX(1,LWORK))
          On exit, if INFO=0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.
          If N <= 1,               LWORK >= 1.
          If JOBZ = 'N' and N > 1, LWORK >= N.
          If JOBZ = 'V' and N > 1, LWORK >= 2*N**2.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal sizes of the WORK, RWORK and
          IWORK arrays, returns these values as the first entries of
          the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]RWORK
          RWORK is REAL array, dimension (MAX(1,LRWORK))
          On exit, if INFO=0, RWORK(1) returns the optimal LRWORK.
[in]LRWORK
          LRWORK is INTEGER
          The dimension of array RWORK.
          If N <= 1,               LRWORK >= 1.
          If JOBZ = 'N' and N > 1, LRWORK >= N.
          If JOBZ = 'V' and N > 1, LRWORK >= 1 + 5*N + 2*N**2.

          If LRWORK = -1, then a workspace query is assumed; the
          routine only calculates the optimal sizes of the WORK, RWORK
          and IWORK arrays, returns these values as the first entries
          of the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO=0, IWORK(1) returns the optimal LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of array IWORK.
          If JOBZ = 'N' or N <= 1, LIWORK >= 1.
          If JOBZ = 'V' and N > 1, LIWORK >= 3 + 5*N.

          If LIWORK = -1, then a workspace query is assumed; the
          routine only calculates the optimal sizes of the WORK, RWORK
          and IWORK arrays, returns these values as the first entries
          of the WORK, RWORK and IWORK arrays, and no error message
          related to LWORK or LRWORK or LIWORK is issued by XERBLA.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  if INFO = i, and i is:
             <= N:  the algorithm failed to converge:
                    i off-diagonal elements of an intermediate
                    tridiagonal form did not converge to zero;
             > N:   if INFO = N + i, for 1 <= i <= N, then CPBSTF
                    returned INFO = i: B is not positive definite.
                    The factorization of B could not be completed and
                    no eigenvalues or eigenvectors were computed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 243 of file chbgvd.f.

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*
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 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
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine cpbstf(uplo, n, kd, ab, ldab, info)
CPBSTF
Definition cpbstf.f:153
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
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
Here is the call graph for this function:
Here is the caller graph for this function: