LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine chbgv ( 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,
real, dimension( * )  RWORK,
integer  INFO 
)

CHBGV

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

Purpose:
 CHBGV 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.
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 (N)
[out]RWORK
          RWORK is REAL array, dimension (3*N)
[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
November 2015

Definition at line 185 of file chbgv.f.

185 *
186 * -- LAPACK driver routine (version 3.6.0) --
187 * -- LAPACK is a software package provided by Univ. of Tennessee, --
188 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
189 * November 2015
190 *
191 * .. Scalar Arguments ..
192  CHARACTER jobz, uplo
193  INTEGER info, ka, kb, ldab, ldbb, ldz, n
194 * ..
195 * .. Array Arguments ..
196  REAL rwork( * ), w( * )
197  COMPLEX ab( ldab, * ), bb( ldbb, * ), work( * ),
198  $ z( ldz, * )
199 * ..
200 *
201 * =====================================================================
202 *
203 * .. Local Scalars ..
204  LOGICAL upper, wantz
205  CHARACTER vect
206  INTEGER iinfo, inde, indwrk
207 * ..
208 * .. External Functions ..
209  LOGICAL lsame
210  EXTERNAL lsame
211 * ..
212 * .. External Subroutines ..
213  EXTERNAL chbgst, chbtrd, cpbstf, csteqr, ssterf, xerbla
214 * ..
215 * .. Executable Statements ..
216 *
217 * Test the input parameters.
218 *
219  wantz = lsame( jobz, 'V' )
220  upper = lsame( uplo, 'U' )
221 *
222  info = 0
223  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
224  info = -1
225  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
226  info = -2
227  ELSE IF( n.LT.0 ) THEN
228  info = -3
229  ELSE IF( ka.LT.0 ) THEN
230  info = -4
231  ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
232  info = -5
233  ELSE IF( ldab.LT.ka+1 ) THEN
234  info = -7
235  ELSE IF( ldbb.LT.kb+1 ) THEN
236  info = -9
237  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
238  info = -12
239  END IF
240  IF( info.NE.0 ) THEN
241  CALL xerbla( 'CHBGV ', -info )
242  RETURN
243  END IF
244 *
245 * Quick return if possible
246 *
247  IF( n.EQ.0 )
248  $ RETURN
249 *
250 * Form a split Cholesky factorization of B.
251 *
252  CALL cpbstf( uplo, n, kb, bb, ldbb, info )
253  IF( info.NE.0 ) THEN
254  info = n + info
255  RETURN
256  END IF
257 *
258 * Transform problem to standard eigenvalue problem.
259 *
260  inde = 1
261  indwrk = inde + n
262  CALL chbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, z, ldz,
263  $ work, rwork( indwrk ), iinfo )
264 *
265 * Reduce to tridiagonal form.
266 *
267  IF( wantz ) THEN
268  vect = 'U'
269  ELSE
270  vect = 'N'
271  END IF
272  CALL chbtrd( vect, uplo, n, ka, ab, ldab, w, rwork( inde ), z,
273  $ ldz, work, iinfo )
274 *
275 * For eigenvalues only, call SSTERF. For eigenvectors, call CSTEQR.
276 *
277  IF( .NOT.wantz ) THEN
278  CALL ssterf( n, w, rwork( inde ), info )
279  ELSE
280  CALL csteqr( jobz, n, w, rwork( inde ), z, ldz,
281  $ rwork( indwrk ), info )
282  END IF
283  RETURN
284 *
285 * End of CHBGV
286 *
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 csteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
CSTEQR
Definition: csteqr.f:134
subroutine chbgst(VECT, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, X, LDX, WORK, RWORK, INFO)
CHBGST
Definition: chbgst.f:167
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:88
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: