LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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.

 The divide and conquer algorithm makes very mild assumptions about
 floating point arithmetic. It will work on machines with a guard
 digit in add/subtract, or on those binary machines without guard
 digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
 Cray-2. It could conceivably fail on hexadecimal or decimal machines
 without guard digits, but we know of none.
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.
Date
June 2016
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 254 of file chbgvd.f.

254 *
255 * -- LAPACK driver routine (version 3.6.1) --
256 * -- LAPACK is a software package provided by Univ. of Tennessee, --
257 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
258 * June 2016
259 *
260 * .. Scalar Arguments ..
261  CHARACTER jobz, uplo
262  INTEGER info, ka, kb, ldab, ldbb, ldz, liwork, lrwork,
263  $ lwork, n
264 * ..
265 * .. Array Arguments ..
266  INTEGER iwork( * )
267  REAL rwork( * ), w( * )
268  COMPLEX ab( ldab, * ), bb( ldbb, * ), work( * ),
269  $ z( ldz, * )
270 * ..
271 *
272 * =====================================================================
273 *
274 * .. Parameters ..
275  COMPLEX cone, czero
276  parameter ( cone = ( 1.0e+0, 0.0e+0 ),
277  $ czero = ( 0.0e+0, 0.0e+0 ) )
278 * ..
279 * .. Local Scalars ..
280  LOGICAL lquery, upper, wantz
281  CHARACTER vect
282  INTEGER iinfo, inde, indwk2, indwrk, liwmin, llrwk,
283  $ llwk2, lrwmin, lwmin
284 * ..
285 * .. External Functions ..
286  LOGICAL lsame
287  EXTERNAL lsame
288 * ..
289 * .. External Subroutines ..
290  EXTERNAL ssterf, xerbla, cgemm, chbgst, chbtrd, clacpy,
291  $ cpbstf, cstedc
292 * ..
293 * .. Executable Statements ..
294 *
295 * Test the input parameters.
296 *
297  wantz = lsame( jobz, 'V' )
298  upper = lsame( uplo, 'U' )
299  lquery = ( lwork.EQ.-1 .OR. lrwork.EQ.-1 .OR. liwork.EQ.-1 )
300 *
301  info = 0
302  IF( n.LE.1 ) THEN
303  lwmin = 1+n
304  lrwmin = 1+n
305  liwmin = 1
306  ELSE IF( wantz ) THEN
307  lwmin = 2*n**2
308  lrwmin = 1 + 5*n + 2*n**2
309  liwmin = 3 + 5*n
310  ELSE
311  lwmin = n
312  lrwmin = n
313  liwmin = 1
314  END IF
315  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
316  info = -1
317  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
318  info = -2
319  ELSE IF( n.LT.0 ) THEN
320  info = -3
321  ELSE IF( ka.LT.0 ) THEN
322  info = -4
323  ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
324  info = -5
325  ELSE IF( ldab.LT.ka+1 ) THEN
326  info = -7
327  ELSE IF( ldbb.LT.kb+1 ) THEN
328  info = -9
329  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
330  info = -12
331  END IF
332 *
333  IF( info.EQ.0 ) THEN
334  work( 1 ) = lwmin
335  rwork( 1 ) = lrwmin
336  iwork( 1 ) = liwmin
337 *
338  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
339  info = -14
340  ELSE IF( lrwork.LT.lrwmin .AND. .NOT.lquery ) THEN
341  info = -16
342  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
343  info = -18
344  END IF
345  END IF
346 *
347  IF( info.NE.0 ) THEN
348  CALL xerbla( 'CHBGVD', -info )
349  RETURN
350  ELSE IF( lquery ) THEN
351  RETURN
352  END IF
353 *
354 * Quick return if possible
355 *
356  IF( n.EQ.0 )
357  $ RETURN
358 *
359 * Form a split Cholesky factorization of B.
360 *
361  CALL cpbstf( uplo, n, kb, bb, ldbb, info )
362  IF( info.NE.0 ) THEN
363  info = n + info
364  RETURN
365  END IF
366 *
367 * Transform problem to standard eigenvalue problem.
368 *
369  inde = 1
370  indwrk = inde + n
371  indwk2 = 1 + n*n
372  llwk2 = lwork - indwk2 + 2
373  llrwk = lrwork - indwrk + 2
374  CALL chbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, z, ldz,
375  $ work, rwork, iinfo )
376 *
377 * Reduce Hermitian band matrix to tridiagonal form.
378 *
379  IF( wantz ) THEN
380  vect = 'U'
381  ELSE
382  vect = 'N'
383  END IF
384  CALL chbtrd( vect, uplo, n, ka, ab, ldab, w, rwork( inde ), z,
385  $ ldz, work, iinfo )
386 *
387 * For eigenvalues only, call SSTERF. For eigenvectors, call CSTEDC.
388 *
389  IF( .NOT.wantz ) THEN
390  CALL ssterf( n, w, rwork( inde ), info )
391  ELSE
392  CALL cstedc( 'I', n, w, rwork( inde ), work, n, work( indwk2 ),
393  $ llwk2, rwork( indwrk ), llrwk, iwork, liwork,
394  $ info )
395  CALL cgemm( 'N', 'N', n, n, n, cone, z, ldz, work, n, czero,
396  $ work( indwk2 ), n )
397  CALL clacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
398  END IF
399 *
400  work( 1 ) = lwmin
401  rwork( 1 ) = lrwmin
402  iwork( 1 ) = liwmin
403  RETURN
404 *
405 * End of CHBGVD
406 *
subroutine cpbstf(UPLO, N, KD, AB, LDAB, INFO)
CPBSTF
Definition: cpbstf.f:155
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine chbtrd(VECT, UPLO, N, KD, AB, LDAB, D, E, Q, LDQ, WORK, INFO)
CHBTRD
Definition: chbtrd.f:165
subroutine chbgst(VECT, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, X, LDX, WORK, RWORK, INFO)
CHBGST
Definition: chbgst.f:167
subroutine cstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, RWORK, LRWORK, IWORK, LIWORK, INFO)
CSTEDC
Definition: cstedc.f:214
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:105
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:88
subroutine cgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
CGEMM
Definition: cgemm.f:189
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: