LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dsygvd ( integer  ITYPE,
character  JOBZ,
character  UPLO,
integer  N,
double precision, dimension( lda, * )  A,
integer  LDA,
double precision, dimension( ldb, * )  B,
integer  LDB,
double precision, dimension( * )  W,
double precision, dimension( * )  WORK,
integer  LWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

DSYGVD

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

Purpose:
 DSYGVD computes all the eigenvalues, and optionally, the eigenvectors
 of a real generalized symmetric-definite eigenproblem, of the form
 A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.  Here A and
 B are assumed to be symmetric 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]ITYPE
          ITYPE is INTEGER
          Specifies the problem type to be solved:
          = 1:  A*x = (lambda)*B*x
          = 2:  A*B*x = (lambda)*x
          = 3:  B*A*x = (lambda)*x
[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,out]A
          A is DOUBLE PRECISION array, dimension (LDA, N)
          On entry, the symmetric matrix A.  If UPLO = 'U', the
          leading N-by-N upper triangular part of A contains the
          upper triangular part of the matrix A.  If UPLO = 'L',
          the leading N-by-N lower triangular part of A contains
          the lower triangular part of the matrix A.

          On exit, if JOBZ = 'V', then if INFO = 0, A contains the
          matrix Z of eigenvectors.  The eigenvectors are normalized
          as follows:
          if ITYPE = 1 or 2, Z**T*B*Z = I;
          if ITYPE = 3, Z**T*inv(B)*Z = I.
          If JOBZ = 'N', then on exit the upper triangle (if UPLO='U')
          or the lower triangle (if UPLO='L') of A, including the
          diagonal, is destroyed.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[in,out]B
          B is DOUBLE PRECISION array, dimension (LDB, N)
          On entry, the symmetric matrix B.  If UPLO = 'U', the
          leading N-by-N upper triangular part of B contains the
          upper triangular part of the matrix B.  If UPLO = 'L',
          the leading N-by-N lower triangular part of B contains
          the lower triangular part of the matrix B.

          On exit, if INFO <= N, the part of B containing the matrix is
          overwritten by the triangular factor U or L from the Cholesky
          factorization B = U**T*U or B = L*L**T.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[out]W
          W is DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[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+1.
          If JOBZ = 'V' and N > 1, LWORK >= 1 + 6*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 INFO = 0, IWORK(1) returns the optimal LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of the array IWORK.
          If N <= 1,                LIWORK >= 1.
          If JOBZ  = 'N' and 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:  DPOTRF or DSYEVD returned an error code:
             <= N:  if INFO = i and JOBZ = 'N', then the algorithm
                    failed to converge; i off-diagonal elements of an
                    intermediate tridiagonal form did not converge to
                    zero;
                    if INFO = i and JOBZ = 'V', then the algorithm
                    failed to compute an eigenvalue while working on
                    the submatrix lying in rows and columns INFO/(N+1)
                    through mod(INFO,N+1);
             > N:   if INFO = N + i, for 1 <= i <= N, then the leading
                    minor of order i of 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
Further Details:
  Modified so that no backsubstitution is performed if DSYEVD fails to
  converge (NEIG in old code could be greater than N causing out of
  bounds reference to A - reported by Ralf Meyer).  Also corrected the
  description of INFO and the test on ITYPE. Sven, 16 Feb 05.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 229 of file dsygvd.f.

229 *
230 * -- LAPACK driver routine (version 3.6.0) --
231 * -- LAPACK is a software package provided by Univ. of Tennessee, --
232 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
233 * November 2015
234 *
235 * .. Scalar Arguments ..
236  CHARACTER jobz, uplo
237  INTEGER info, itype, lda, ldb, liwork, lwork, n
238 * ..
239 * .. Array Arguments ..
240  INTEGER iwork( * )
241  DOUBLE PRECISION a( lda, * ), b( ldb, * ), w( * ), work( * )
242 * ..
243 *
244 * =====================================================================
245 *
246 * .. Parameters ..
247  DOUBLE PRECISION one
248  parameter ( one = 1.0d+0 )
249 * ..
250 * .. Local Scalars ..
251  LOGICAL lquery, upper, wantz
252  CHARACTER trans
253  INTEGER liopt, liwmin, lopt, lwmin
254 * ..
255 * .. External Functions ..
256  LOGICAL lsame
257  EXTERNAL lsame
258 * ..
259 * .. External Subroutines ..
260  EXTERNAL dpotrf, dsyevd, dsygst, dtrmm, dtrsm, xerbla
261 * ..
262 * .. Intrinsic Functions ..
263  INTRINSIC dble, max
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 + 6*n + 2*n**2
280  ELSE
281  liwmin = 1
282  lwmin = 2*n + 1
283  END IF
284  lopt = lwmin
285  liopt = liwmin
286  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
287  info = -1
288  ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
289  info = -2
290  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
291  info = -3
292  ELSE IF( n.LT.0 ) THEN
293  info = -4
294  ELSE IF( lda.LT.max( 1, n ) ) THEN
295  info = -6
296  ELSE IF( ldb.LT.max( 1, n ) ) THEN
297  info = -8
298  END IF
299 *
300  IF( info.EQ.0 ) THEN
301  work( 1 ) = lopt
302  iwork( 1 ) = liopt
303 *
304  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
305  info = -11
306  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
307  info = -13
308  END IF
309  END IF
310 *
311  IF( info.NE.0 ) THEN
312  CALL xerbla( 'DSYGVD', -info )
313  RETURN
314  ELSE IF( lquery ) THEN
315  RETURN
316  END IF
317 *
318 * Quick return if possible
319 *
320  IF( n.EQ.0 )
321  $ RETURN
322 *
323 * Form a Cholesky factorization of B.
324 *
325  CALL dpotrf( uplo, n, b, ldb, info )
326  IF( info.NE.0 ) THEN
327  info = n + info
328  RETURN
329  END IF
330 *
331 * Transform problem to standard eigenvalue problem and solve.
332 *
333  CALL dsygst( itype, uplo, n, a, lda, b, ldb, info )
334  CALL dsyevd( jobz, uplo, n, a, lda, w, work, lwork, iwork, liwork,
335  $ info )
336  lopt = max( dble( lopt ), dble( work( 1 ) ) )
337  liopt = max( dble( liopt ), dble( iwork( 1 ) ) )
338 *
339  IF( wantz .AND. info.EQ.0 ) THEN
340 *
341 * Backtransform eigenvectors to the original problem.
342 *
343  IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
344 *
345 * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
346 * backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
347 *
348  IF( upper ) THEN
349  trans = 'N'
350  ELSE
351  trans = 'T'
352  END IF
353 *
354  CALL dtrsm( 'Left', uplo, trans, 'Non-unit', n, n, one,
355  $ b, ldb, a, lda )
356 *
357  ELSE IF( itype.EQ.3 ) THEN
358 *
359 * For B*A*x=(lambda)*x;
360 * backtransform eigenvectors: x = L*y or U**T*y
361 *
362  IF( upper ) THEN
363  trans = 'T'
364  ELSE
365  trans = 'N'
366  END IF
367 *
368  CALL dtrmm( 'Left', uplo, trans, 'Non-unit', n, n, one,
369  $ b, ldb, a, lda )
370  END IF
371  END IF
372 *
373  work( 1 ) = lopt
374  iwork( 1 ) = liopt
375 *
376  RETURN
377 *
378 * End of DSYGVD
379 *
subroutine dtrsm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
DTRSM
Definition: dtrsm.f:183
subroutine dpotrf(UPLO, N, A, LDA, INFO)
DPOTRF
Definition: dpotrf.f:109
subroutine dtrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
DTRMM
Definition: dtrmm.f:179
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dsyevd(JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, IWORK, LIWORK, INFO)
DSYEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices ...
Definition: dsyevd.f:187
subroutine dsygst(ITYPE, UPLO, N, A, LDA, B, LDB, INFO)
DSYGST
Definition: dsygst.f:129
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: