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

◆ dsygvd()

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.
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
                    principal minor of order i of B is not positive.
                    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.
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 219 of file dsygvd.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, ITYPE, LDA, LDB, LIWORK, LWORK, N
229* ..
230* .. Array Arguments ..
231 INTEGER IWORK( * )
232 DOUBLE PRECISION A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
233* ..
234*
235* =====================================================================
236*
237* .. Parameters ..
238 DOUBLE PRECISION ONE
239 parameter( one = 1.0d+0 )
240* ..
241* .. Local Scalars ..
242 LOGICAL LQUERY, UPPER, WANTZ
243 CHARACTER TRANS
244 INTEGER LIOPT, LIWMIN, LOPT, LWMIN
245* ..
246* .. External Functions ..
247 LOGICAL LSAME
248 EXTERNAL lsame
249* ..
250* .. External Subroutines ..
251 EXTERNAL dpotrf, dsyevd, dsygst, dtrmm, dtrsm, xerbla
252* ..
253* .. Intrinsic Functions ..
254 INTRINSIC dble, max
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 + 6*n + 2*n**2
271 ELSE
272 liwmin = 1
273 lwmin = 2*n + 1
274 END IF
275 lopt = lwmin
276 liopt = liwmin
277 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
278 info = -1
279 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
280 info = -2
281 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
282 info = -3
283 ELSE IF( n.LT.0 ) THEN
284 info = -4
285 ELSE IF( lda.LT.max( 1, n ) ) THEN
286 info = -6
287 ELSE IF( ldb.LT.max( 1, n ) ) THEN
288 info = -8
289 END IF
290*
291 IF( info.EQ.0 ) THEN
292 work( 1 ) = lopt
293 iwork( 1 ) = liopt
294*
295 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
296 info = -11
297 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
298 info = -13
299 END IF
300 END IF
301*
302 IF( info.NE.0 ) THEN
303 CALL xerbla( 'DSYGVD', -info )
304 RETURN
305 ELSE IF( lquery ) THEN
306 RETURN
307 END IF
308*
309* Quick return if possible
310*
311 IF( n.EQ.0 )
312 $ RETURN
313*
314* Form a Cholesky factorization of B.
315*
316 CALL dpotrf( uplo, n, b, ldb, info )
317 IF( info.NE.0 ) THEN
318 info = n + info
319 RETURN
320 END IF
321*
322* Transform problem to standard eigenvalue problem and solve.
323*
324 CALL dsygst( itype, uplo, n, a, lda, b, ldb, info )
325 CALL dsyevd( jobz, uplo, n, a, lda, w, work, lwork, iwork, liwork,
326 $ info )
327 lopt = int( max( dble( lopt ), dble( work( 1 ) ) ) )
328 liopt = int( max( dble( liopt ), dble( iwork( 1 ) ) ) )
329*
330 IF( wantz .AND. info.EQ.0 ) THEN
331*
332* Backtransform eigenvectors to the original problem.
333*
334 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
335*
336* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
337* backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
338*
339 IF( upper ) THEN
340 trans = 'N'
341 ELSE
342 trans = 'T'
343 END IF
344*
345 CALL dtrsm( 'Left', uplo, trans, 'Non-unit', n, n, one,
346 $ b, ldb, a, lda )
347*
348 ELSE IF( itype.EQ.3 ) THEN
349*
350* For B*A*x=(lambda)*x;
351* backtransform eigenvectors: x = L*y or U**T*y
352*
353 IF( upper ) THEN
354 trans = 'T'
355 ELSE
356 trans = 'N'
357 END IF
358*
359 CALL dtrmm( 'Left', uplo, trans, 'Non-unit', n, n, one,
360 $ b, ldb, a, lda )
361 END IF
362 END IF
363*
364 work( 1 ) = lopt
365 iwork( 1 ) = liopt
366*
367 RETURN
368*
369* End of DSYGVD
370*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
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:178
subroutine dsygst(itype, uplo, n, a, lda, b, ldb, info)
DSYGST
Definition dsygst.f:127
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dpotrf(uplo, n, a, lda, info)
DPOTRF
Definition dpotrf.f:107
subroutine dtrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRMM
Definition dtrmm.f:177
subroutine dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
DTRSM
Definition dtrsm.f:181
Here is the call graph for this function:
Here is the caller graph for this function: