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

DSYEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matrices

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

Purpose:
 DSYEVD computes all eigenvalues and, optionally, eigenvectors of a
 real symmetric matrix A. 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.

 Because of large use of BLAS of level 3, DSYEVD needs N**2 more
 workspace than DSYEVX.
Parameters
[in]JOBZ
          JOBZ is CHARACTER*1
          = 'N':  Compute eigenvalues only;
          = 'V':  Compute eigenvalues and eigenvectors.
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  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
          orthonormal eigenvectors of the matrix A.
          If JOBZ = 'N', then on exit the lower triangle (if UPLO='L')
          or the upper triangle (if UPLO='U') of A, including the
          diagonal, is destroyed.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= 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 (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 must be at least 1.
          If JOBZ = 'N' and N > 1, LWORK must be at least 2*N+1.
          If JOBZ = 'V' and N > 1, LWORK must be at least
                                                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 must be at least 1.
          If JOBZ  = 'N' and N > 1, LIWORK must be at least 1.
          If JOBZ  = 'V' and N > 1, LIWORK must be at least 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 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).
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012
Contributors:
Jeff Rutter, Computer Science Division, University of California at Berkeley, USA
Modified by Francoise Tisseur, University of Tennessee
Modified description of INFO. Sven, 16 Feb 05.

Definition at line 187 of file dsyevd.f.

187 *
188 * -- LAPACK driver routine (version 3.4.2) --
189 * -- LAPACK is a software package provided by Univ. of Tennessee, --
190 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
191 * September 2012
192 *
193 * .. Scalar Arguments ..
194  CHARACTER jobz, uplo
195  INTEGER info, lda, liwork, lwork, n
196 * ..
197 * .. Array Arguments ..
198  INTEGER iwork( * )
199  DOUBLE PRECISION a( lda, * ), w( * ), work( * )
200 * ..
201 *
202 * =====================================================================
203 *
204 * .. Parameters ..
205  DOUBLE PRECISION zero, one
206  parameter ( zero = 0.0d+0, one = 1.0d+0 )
207 * ..
208 * .. Local Scalars ..
209 *
210  LOGICAL lower, lquery, wantz
211  INTEGER iinfo, inde, indtau, indwk2, indwrk, iscale,
212  $ liopt, liwmin, llwork, llwrk2, lopt, lwmin
213  DOUBLE PRECISION anrm, bignum, eps, rmax, rmin, safmin, sigma,
214  $ smlnum
215 * ..
216 * .. External Functions ..
217  LOGICAL lsame
218  INTEGER ilaenv
219  DOUBLE PRECISION dlamch, dlansy
220  EXTERNAL lsame, dlamch, dlansy, ilaenv
221 * ..
222 * .. External Subroutines ..
223  EXTERNAL dlacpy, dlascl, dormtr, dscal, dstedc, dsterf,
224  $ dsytrd, xerbla
225 * ..
226 * .. Intrinsic Functions ..
227  INTRINSIC max, sqrt
228 * ..
229 * .. Executable Statements ..
230 *
231 * Test the input parameters.
232 *
233  wantz = lsame( jobz, 'V' )
234  lower = lsame( uplo, 'L' )
235  lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
236 *
237  info = 0
238  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
239  info = -1
240  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
241  info = -2
242  ELSE IF( n.LT.0 ) THEN
243  info = -3
244  ELSE IF( lda.LT.max( 1, n ) ) THEN
245  info = -5
246  END IF
247 *
248  IF( info.EQ.0 ) THEN
249  IF( n.LE.1 ) THEN
250  liwmin = 1
251  lwmin = 1
252  lopt = lwmin
253  liopt = liwmin
254  ELSE
255  IF( wantz ) THEN
256  liwmin = 3 + 5*n
257  lwmin = 1 + 6*n + 2*n**2
258  ELSE
259  liwmin = 1
260  lwmin = 2*n + 1
261  END IF
262  lopt = max( lwmin, 2*n +
263  $ ilaenv( 1, 'DSYTRD', uplo, n, -1, -1, -1 ) )
264  liopt = liwmin
265  END IF
266  work( 1 ) = lopt
267  iwork( 1 ) = liopt
268 *
269  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
270  info = -8
271  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
272  info = -10
273  END IF
274  END IF
275 *
276  IF( info.NE.0 ) THEN
277  CALL xerbla( 'DSYEVD', -info )
278  RETURN
279  ELSE IF( lquery ) THEN
280  RETURN
281  END IF
282 *
283 * Quick return if possible
284 *
285  IF( n.EQ.0 )
286  $ RETURN
287 *
288  IF( n.EQ.1 ) THEN
289  w( 1 ) = a( 1, 1 )
290  IF( wantz )
291  $ a( 1, 1 ) = one
292  RETURN
293  END IF
294 *
295 * Get machine constants.
296 *
297  safmin = dlamch( 'Safe minimum' )
298  eps = dlamch( 'Precision' )
299  smlnum = safmin / eps
300  bignum = one / smlnum
301  rmin = sqrt( smlnum )
302  rmax = sqrt( bignum )
303 *
304 * Scale matrix to allowable range, if necessary.
305 *
306  anrm = dlansy( 'M', uplo, n, a, lda, work )
307  iscale = 0
308  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
309  iscale = 1
310  sigma = rmin / anrm
311  ELSE IF( anrm.GT.rmax ) THEN
312  iscale = 1
313  sigma = rmax / anrm
314  END IF
315  IF( iscale.EQ.1 )
316  $ CALL dlascl( uplo, 0, 0, one, sigma, n, n, a, lda, info )
317 *
318 * Call DSYTRD to reduce symmetric matrix to tridiagonal form.
319 *
320  inde = 1
321  indtau = inde + n
322  indwrk = indtau + n
323  llwork = lwork - indwrk + 1
324  indwk2 = indwrk + n*n
325  llwrk2 = lwork - indwk2 + 1
326 *
327  CALL dsytrd( uplo, n, a, lda, w, work( inde ), work( indtau ),
328  $ work( indwrk ), llwork, iinfo )
329 *
330 * For eigenvalues only, call DSTERF. For eigenvectors, first call
331 * DSTEDC to generate the eigenvector matrix, WORK(INDWRK), of the
332 * tridiagonal matrix, then call DORMTR to multiply it by the
333 * Householder transformations stored in A.
334 *
335  IF( .NOT.wantz ) THEN
336  CALL dsterf( n, w, work( inde ), info )
337  ELSE
338  CALL dstedc( 'I', n, w, work( inde ), work( indwrk ), n,
339  $ work( indwk2 ), llwrk2, iwork, liwork, info )
340  CALL dormtr( 'L', uplo, 'N', n, n, a, lda, work( indtau ),
341  $ work( indwrk ), n, work( indwk2 ), llwrk2, iinfo )
342  CALL dlacpy( 'A', n, n, work( indwrk ), n, a, lda )
343  END IF
344 *
345 * If matrix was scaled, then rescale eigenvalues appropriately.
346 *
347  IF( iscale.EQ.1 )
348  $ CALL dscal( n, one / sigma, w, 1 )
349 *
350  work( 1 ) = lopt
351  iwork( 1 ) = liopt
352 *
353  RETURN
354 *
355 * End of DSYEVD
356 *
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
double precision function dlansy(NORM, UPLO, N, A, LDA, WORK)
DLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a real symmetric matrix.
Definition: dlansy.f:124
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dsytrd(UPLO, N, A, LDA, D, E, TAU, WORK, LWORK, INFO)
DSYTRD
Definition: dsytrd.f:194
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dlascl(TYPE, KL, KU, CFROM, CTO, M, N, A, LDA, INFO)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition: dlascl.f:145
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dormtr(SIDE, UPLO, TRANS, M, N, A, LDA, TAU, C, LDC, WORK, LWORK, INFO)
DORMTR
Definition: dormtr.f:173
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
subroutine dstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
DSTEDC
Definition: dstedc.f:191
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
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: