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

◆ dsbgvd()

subroutine dsbgvd ( 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  lwork,
integer, dimension( * )  iwork,
integer  liwork,
integer  info 
)

DSBGVD

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

Purpose:
 DSBGVD 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 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(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 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 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 DOUBLE PRECISION 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 >= 2*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 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.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 219 of file dsbgvd.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 DOUBLE PRECISION AB( LDAB, * ), BB( LDBB, * ), W( * ),
233 $ WORK( * ), Z( LDZ, * )
234* ..
235*
236* =====================================================================
237*
238* .. Parameters ..
239 DOUBLE PRECISION ONE, ZERO
240 parameter( one = 1.0d+0, zero = 0.0d+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 EXTERNAL lsame
251* ..
252* .. External Subroutines ..
253 EXTERNAL dgemm, dlacpy, dpbstf, dsbgst, dsbtrd, dstedc,
254 $ dsterf, xerbla
255* ..
256* .. Executable Statements ..
257*
258* Test the input parameters.
259*
260 wantz = lsame( jobz, 'V' )
261 upper = lsame( uplo, 'U' )
262 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
263*
264 info = 0
265 IF( n.LE.1 ) THEN
266 liwmin = 1
267 lwmin = 1
268 ELSE IF( wantz ) THEN
269 liwmin = 3 + 5*n
270 lwmin = 1 + 5*n + 2*n**2
271 ELSE
272 liwmin = 1
273 lwmin = 2*n
274 END IF
275*
276 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
277 info = -1
278 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
279 info = -2
280 ELSE IF( n.LT.0 ) THEN
281 info = -3
282 ELSE IF( ka.LT.0 ) THEN
283 info = -4
284 ELSE IF( kb.LT.0 .OR. kb.GT.ka ) THEN
285 info = -5
286 ELSE IF( ldab.LT.ka+1 ) THEN
287 info = -7
288 ELSE IF( ldbb.LT.kb+1 ) THEN
289 info = -9
290 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
291 info = -12
292 END IF
293*
294 IF( info.EQ.0 ) THEN
295 work( 1 ) = lwmin
296 iwork( 1 ) = liwmin
297*
298 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
299 info = -14
300 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
301 info = -16
302 END IF
303 END IF
304*
305 IF( info.NE.0 ) THEN
306 CALL xerbla( 'DSBGVD', -info )
307 RETURN
308 ELSE IF( lquery ) THEN
309 RETURN
310 END IF
311*
312* Quick return if possible
313*
314 IF( n.EQ.0 )
315 $ RETURN
316*
317* Form a split Cholesky factorization of B.
318*
319 CALL dpbstf( uplo, n, kb, bb, ldbb, info )
320 IF( info.NE.0 ) THEN
321 info = n + info
322 RETURN
323 END IF
324*
325* Transform problem to standard eigenvalue problem.
326*
327 inde = 1
328 indwrk = inde + n
329 indwk2 = indwrk + n*n
330 llwrk2 = lwork - indwk2 + 1
331 CALL dsbgst( jobz, uplo, n, ka, kb, ab, ldab, bb, ldbb, z, ldz,
332 $ work, iinfo )
333*
334* Reduce to tridiagonal form.
335*
336 IF( wantz ) THEN
337 vect = 'U'
338 ELSE
339 vect = 'N'
340 END IF
341 CALL dsbtrd( vect, uplo, n, ka, ab, ldab, w, work( inde ), z, ldz,
342 $ work( indwrk ), iinfo )
343*
344* For eigenvalues only, call DSTERF. For eigenvectors, call SSTEDC.
345*
346 IF( .NOT.wantz ) THEN
347 CALL dsterf( n, w, work( inde ), info )
348 ELSE
349 CALL dstedc( 'I', n, w, work( inde ), work( indwrk ), n,
350 $ work( indwk2 ), llwrk2, iwork, liwork, info )
351 CALL dgemm( 'N', 'N', n, n, n, one, z, ldz, work( indwrk ), n,
352 $ zero, work( indwk2 ), n )
353 CALL dlacpy( 'A', n, n, work( indwk2 ), n, z, ldz )
354 END IF
355*
356 work( 1 ) = lwmin
357 iwork( 1 ) = liwmin
358*
359 RETURN
360*
361* End of DSBGVD
362*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
DGEMM
Definition dgemm.f:188
subroutine dsbgst(vect, uplo, n, ka, kb, ab, ldab, bb, ldbb, x, ldx, work, info)
DSBGST
Definition dsbgst.f:159
subroutine dsbtrd(vect, uplo, n, kd, ab, ldab, d, e, q, ldq, work, info)
DSBTRD
Definition dsbtrd.f:163
subroutine dlacpy(uplo, m, n, a, lda, b, ldb)
DLACPY copies all or part of one two-dimensional array to another.
Definition dlacpy.f:103
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dpbstf(uplo, n, kd, ab, ldab, info)
DPBSTF
Definition dpbstf.f:152
subroutine dstedc(compz, n, d, e, z, ldz, work, lwork, iwork, liwork, info)
DSTEDC
Definition dstedc.f:182
subroutine dsterf(n, d, e, info)
DSTERF
Definition dsterf.f:86
Here is the call graph for this function:
Here is the caller graph for this function: