LAPACK  3.10.1
LAPACK: Linear Algebra PACKage

◆ chbgv()

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.

Definition at line 181 of file chbgv.f.

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