LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine ssbgvd ( character  JOBZ,
character  UPLO,
integer  N,
integer  KA,
integer  KB,
real, dimension( ldab, * )  AB,
integer  LDAB,
real, dimension( ldbb, * )  BB,
integer  LDBB,
real, dimension( * )  W,
real, dimension( ldz, * )  Z,
integer  LDZ,
real, dimension( * )  WORK,
integer  LWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

SSBGVD

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

Purpose:
 SSBGVD 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.  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 REAL 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 REAL 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(ka+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 SPBSTF.
[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 REAL 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 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 >= max(1,N).
[out]WORK
          WORK is REAL 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 >= 3*N.
          If JOBZ = 'V' and N > 1, LWORK >= 1 + 5*N + 2*N**2.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal sizes of the WORK and IWORK
          arrays, returns these values as the first entries of the WORK
          and IWORK arrays, and no error message related to LWORK or
          LIWORK is issued by XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if LIWORK > 0, IWORK(1) returns the optimal LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of the 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 and
          IWORK arrays, returns these values as the first entries of
          the WORK and IWORK arrays, and no error message related to
          LWORK 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 SPBSTF
                    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 229 of file ssbgvd.f.

229 *
230 * -- LAPACK driver routine (version 3.6.1) --
231 * -- LAPACK is a software package provided by Univ. of Tennessee, --
232 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
233 * June 2016
234 *
235 * .. Scalar Arguments ..
236  CHARACTER jobz, uplo
237  INTEGER info, ka, kb, ldab, ldbb, ldz, liwork, lwork, n
238 * ..
239 * .. Array Arguments ..
240  INTEGER iwork( * )
241  REAL ab( ldab, * ), bb( ldbb, * ), w( * ),
242  $ work( * ), z( ldz, * )
243 * ..
244 *
245 * =====================================================================
246 *
247 * .. Parameters ..
248  REAL one, zero
249  parameter ( one = 1.0e+0, zero = 0.0e+0 )
250 * ..
251 * .. Local Scalars ..
252  LOGICAL lquery, upper, wantz
253  CHARACTER vect
254  INTEGER iinfo, inde, indwk2, indwrk, liwmin, llwrk2,
255  $ lwmin
256 * ..
257 * .. External Functions ..
258  LOGICAL lsame
259  EXTERNAL lsame
260 * ..
261 * .. External Subroutines ..
262  EXTERNAL sgemm, slacpy, spbstf, ssbgst, ssbtrd, sstedc,
263  $ ssterf, xerbla
264 * ..
265 * .. Executable Statements ..
266 *
267 * Test the input parameters.
268 *
269  wantz = lsame( jobz, 'V' )
270  upper = lsame( uplo, 'U' )
271  lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
272 *
273  info = 0
274  IF( n.LE.1 ) THEN
275  liwmin = 1
276  lwmin = 1
277  ELSE IF( wantz ) THEN
278  liwmin = 3 + 5*n
279  lwmin = 1 + 5*n + 2*n**2
280  ELSE
281  liwmin = 1
282  lwmin = 2*n
283  END IF
284 *
285  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
286  info = -1
287  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
288  info = -2
289  ELSE IF( n.LT.0 ) THEN
290  info = -3
291  ELSE IF( ka.LT.0 ) THEN
292  info = -4
293  ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
294  info = -5
295  ELSE IF( ldab.LT.ka+1 ) THEN
296  info = -7
297  ELSE IF( ldbb.LT.kb+1 ) THEN
298  info = -9
299  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
300  info = -12
301  END IF
302 *
303  IF( info.EQ.0 ) THEN
304  work( 1 ) = lwmin
305  iwork( 1 ) = liwmin
306 *
307  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
308  info = -14
309  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
310  info = -16
311  END IF
312  END IF
313 *
314  IF( info.NE.0 ) THEN
315  CALL xerbla( 'SSBGVD', -info )
316  RETURN
317  ELSE IF( lquery ) THEN
318  RETURN
319  END IF
320 *
321 * Quick return if possible
322 *
323  IF( n.EQ.0 )
324  $ RETURN
325 *
326 * Form a split Cholesky factorization of B.
327 *
328  CALL spbstf( uplo, n, kb, bb, ldbb, info )
329  IF( info.NE.0 ) THEN
330  info = n + info
331  RETURN
332  END IF
333 *
334 * Transform problem to standard eigenvalue problem.
335 *
336  inde = 1
337  indwrk = inde + n
338  indwk2 = indwrk + n*n
339  llwrk2 = lwork - indwk2 + 1
340  CALL ssbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, z, ldz,
341  $ work, iinfo )
342 *
343 * Reduce to tridiagonal form.
344 *
345  IF( wantz ) THEN
346  vect = 'U'
347  ELSE
348  vect = 'N'
349  END IF
350  CALL ssbtrd( vect, uplo, n, ka, ab, ldab, w, work( inde ), z, ldz,
351  $ work( indwrk ), iinfo )
352 *
353 * For eigenvalues only, call SSTERF. For eigenvectors, call SSTEDC.
354 *
355  IF( .NOT.wantz ) THEN
356  CALL ssterf( n, w, work( inde ), info )
357  ELSE
358  CALL sstedc( 'I', n, w, work( inde ), work( indwrk ), n,
359  $ work( indwk2 ), llwrk2, iwork, liwork, info )
360  CALL sgemm( 'N', 'N', n, n, n, one, z, ldz, work( indwrk ), n,
361  $ zero, work( indwk2 ), n )
362  CALL slacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
363  END IF
364 *
365  work( 1 ) = lwmin
366  iwork( 1 ) = liwmin
367 *
368  RETURN
369 *
370 * End of SSBGVD
371 *
subroutine sgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
SGEMM
Definition: sgemm.f:189
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
subroutine spbstf(UPLO, N, KD, AB, LDAB, INFO)
SPBSTF
Definition: spbstf.f:154
subroutine ssbtrd(VECT, UPLO, N, KD, AB, LDAB, D, E, Q, LDQ, WORK, INFO)
SSBTRD
Definition: ssbtrd.f:165
subroutine ssterf(N, D, E, INFO)
SSTERF
Definition: ssterf.f:88
subroutine sstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
SSTEDC
Definition: sstedc.f:190
subroutine ssbgst(VECT, UPLO, N, KA, KB, AB, LDAB, BB, LDBB, X, LDX, WORK, INFO)
SSBGST
Definition: ssbgst.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: