LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dsbgv ( character  JOBZ,
character  UPLO,
integer  N,
integer  KA,
integer  KB,
double precision, dimension( ldab, * )  AB,
integer  LDAB,
double precision, dimension( ldbb, * )  BB,
integer  LDBB,
double precision, dimension( * )  W,
double precision, dimension( ldz, * )  Z,
integer  LDZ,
double precision, dimension( * )  WORK,
integer  INFO 
)

DSBGV

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

Purpose:
 DSBGV computes all the eigenvalues, and optionally, the eigenvectors
 of a real generalized symmetric-definite banded eigenproblem, of
 the form A*x=(lambda)*B*x. Here A and B are assumed to be symmetric
 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 DOUBLE PRECISION array, dimension (LDAB, N)
          On entry, the upper or lower triangle of the symmetric 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 DOUBLE PRECISION array, dimension (LDBB, N)
          On entry, the upper or lower triangle of the symmetric 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**T*S, as returned by DPBSTF.
[in]LDBB
          LDBB is INTEGER
          The leading dimension of the array BB.  LDBB >= KB+1.
[out]W
          W is DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is DOUBLE PRECISION 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**T*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 DOUBLE PRECISION 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 DPBSTF
                    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 179 of file dsbgv.f.

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