LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ ssbgvd()

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.
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.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 219 of file ssbgvd.f.

221*
222* -- LAPACK driver routine --
223* -- LAPACK is a software package provided by Univ. of Tennessee, --
224* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
225*
226* .. Scalar Arguments ..
227 CHARACTER JOBZ, UPLO
228 INTEGER INFO, KA, KB, LDAB, LDBB, LDZ, LIWORK, LWORK, N
229* ..
230* .. Array Arguments ..
231 INTEGER IWORK( * )
232 REAL AB( LDAB, * ), BB( LDBB, * ), W( * ),
233 $ WORK( * ), Z( LDZ, * )
234* ..
235*
236* =====================================================================
237*
238* .. Parameters ..
239 REAL ONE, ZERO
240 parameter( one = 1.0e+0, zero = 0.0e+0 )
241* ..
242* .. Local Scalars ..
243 LOGICAL LQUERY, UPPER, WANTZ
244 CHARACTER VECT
245 INTEGER IINFO, INDE, INDWK2, INDWRK, LIWMIN, LLWRK2,
246 $ LWMIN
247* ..
248* .. External Functions ..
249 LOGICAL LSAME
250 REAL SROUNDUP_LWORK
251 EXTERNAL lsame, sroundup_lwork
252* ..
253* .. External Subroutines ..
254 EXTERNAL sgemm, slacpy, spbstf, ssbgst, ssbtrd, sstedc,
255 $ ssterf, xerbla
256* ..
257* .. Executable Statements ..
258*
259* Test the input parameters.
260*
261 wantz = lsame( jobz, 'V' )
262 upper = lsame( uplo, 'U' )
263 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
264*
265 info = 0
266 IF( n.LE.1 ) THEN
267 liwmin = 1
268 lwmin = 1
269 ELSE IF( wantz ) THEN
270 liwmin = 3 + 5*n
271 lwmin = 1 + 5*n + 2*n**2
272 ELSE
273 liwmin = 1
274 lwmin = 2*n
275 END IF
276*
277 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
278 info = -1
279 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
280 info = -2
281 ELSE IF( n.LT.0 ) THEN
282 info = -3
283 ELSE IF( ka.LT.0 ) THEN
284 info = -4
285 ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
286 info = -5
287 ELSE IF( ldab.LT.ka+1 ) THEN
288 info = -7
289 ELSE IF( ldbb.LT.kb+1 ) THEN
290 info = -9
291 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
292 info = -12
293 END IF
294*
295 IF( info.EQ.0 ) THEN
296 work( 1 ) = sroundup_lwork(lwmin)
297 iwork( 1 ) = liwmin
298*
299 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
300 info = -14
301 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
302 info = -16
303 END IF
304 END IF
305*
306 IF( info.NE.0 ) THEN
307 CALL xerbla( 'SSBGVD', -info )
308 RETURN
309 ELSE IF( lquery ) THEN
310 RETURN
311 END IF
312*
313* Quick return if possible
314*
315 IF( n.EQ.0 )
316 $ RETURN
317*
318* Form a split Cholesky factorization of B.
319*
320 CALL spbstf( uplo, n, kb, bb, ldbb, info )
321 IF( info.NE.0 ) THEN
322 info = n + info
323 RETURN
324 END IF
325*
326* Transform problem to standard eigenvalue problem.
327*
328 inde = 1
329 indwrk = inde + n
330 indwk2 = indwrk + n*n
331 llwrk2 = lwork - indwk2 + 1
332 CALL ssbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, z, ldz,
333 $ work, iinfo )
334*
335* Reduce to tridiagonal form.
336*
337 IF( wantz ) THEN
338 vect = 'U'
339 ELSE
340 vect = 'N'
341 END IF
342 CALL ssbtrd( vect, uplo, n, ka, ab, ldab, w, work( inde ), z, ldz,
343 $ work( indwrk ), iinfo )
344*
345* For eigenvalues only, call SSTERF. For eigenvectors, call SSTEDC.
346*
347 IF( .NOT.wantz ) THEN
348 CALL ssterf( n, w, work( inde ), info )
349 ELSE
350 CALL sstedc( 'I', n, w, work( inde ), work( indwrk ), n,
351 $ work( indwk2 ), llwrk2, iwork, liwork, info )
352 CALL sgemm( 'N', 'N', n, n, n, one, z, ldz, work( indwrk ), n,
353 $ zero, work( indwk2 ), n )
354 CALL slacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
355 END IF
356*
357 work( 1 ) = sroundup_lwork(lwmin)
358 iwork( 1 ) = liwmin
359*
360 RETURN
361*
362* End of SSBGVD
363*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine sgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SGEMM
Definition sgemm.f:188
subroutine ssbgst(vect, uplo, n, ka, kb, ab, ldab, bb, ldbb, x, ldx, work, info)
SSBGST
Definition ssbgst.f:159
subroutine ssbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
SSBTRD
Definition ssbtrd.f:163
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine spbstf(uplo, n, kd, ab, ldab, info)
SPBSTF
Definition spbstf.f:152
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine sstedc(compz, n, d, e, z, ldz, work, lwork, iwork, liwork, info)
SSTEDC
Definition sstedc.f:182
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:86
Here is the call graph for this function:
Here is the caller graph for this function: