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

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

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

Purpose:
 DSYEV computes all eigenvalues and, optionally, eigenvectors of a
 real symmetric matrix A.
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 (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The length of the array WORK.  LWORK >= max(1,3*N-1).
          For optimal efficiency, LWORK >= (NB+2)*N,
          where NB is the blocksize for DSYTRD returned by ILAENV.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK 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, the algorithm failed to converge; i
                off-diagonal elements of an intermediate tridiagonal
                form did not converge to zero.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 134 of file dsyev.f.

134 *
135 * -- LAPACK driver routine (version 3.4.0) --
136 * -- LAPACK is a software package provided by Univ. of Tennessee, --
137 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
138 * November 2011
139 *
140 * .. Scalar Arguments ..
141  CHARACTER jobz, uplo
142  INTEGER info, lda, lwork, n
143 * ..
144 * .. Array Arguments ..
145  DOUBLE PRECISION a( lda, * ), w( * ), work( * )
146 * ..
147 *
148 * =====================================================================
149 *
150 * .. Parameters ..
151  DOUBLE PRECISION zero, one
152  parameter ( zero = 0.0d0, one = 1.0d0 )
153 * ..
154 * .. Local Scalars ..
155  LOGICAL lower, lquery, wantz
156  INTEGER iinfo, imax, inde, indtau, indwrk, iscale,
157  $ llwork, lwkopt, nb
158  DOUBLE PRECISION anrm, bignum, eps, rmax, rmin, safmin, sigma,
159  $ smlnum
160 * ..
161 * .. External Functions ..
162  LOGICAL lsame
163  INTEGER ilaenv
164  DOUBLE PRECISION dlamch, dlansy
165  EXTERNAL lsame, ilaenv, dlamch, dlansy
166 * ..
167 * .. External Subroutines ..
168  EXTERNAL dlascl, dorgtr, dscal, dsteqr, dsterf, dsytrd,
169  $ xerbla
170 * ..
171 * .. Intrinsic Functions ..
172  INTRINSIC max, sqrt
173 * ..
174 * .. Executable Statements ..
175 *
176 * Test the input parameters.
177 *
178  wantz = lsame( jobz, 'V' )
179  lower = lsame( uplo, 'L' )
180  lquery = ( lwork.EQ.-1 )
181 *
182  info = 0
183  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
184  info = -1
185  ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
186  info = -2
187  ELSE IF( n.LT.0 ) THEN
188  info = -3
189  ELSE IF( lda.LT.max( 1, n ) ) THEN
190  info = -5
191  END IF
192 *
193  IF( info.EQ.0 ) THEN
194  nb = ilaenv( 1, 'DSYTRD', uplo, n, -1, -1, -1 )
195  lwkopt = max( 1, ( nb+2 )*n )
196  work( 1 ) = lwkopt
197 *
198  IF( lwork.LT.max( 1, 3*n-1 ) .AND. .NOT.lquery )
199  $ info = -8
200  END IF
201 *
202  IF( info.NE.0 ) THEN
203  CALL xerbla( 'DSYEV ', -info )
204  RETURN
205  ELSE IF( lquery ) THEN
206  RETURN
207  END IF
208 *
209 * Quick return if possible
210 *
211  IF( n.EQ.0 ) THEN
212  RETURN
213  END IF
214 *
215  IF( n.EQ.1 ) THEN
216  w( 1 ) = a( 1, 1 )
217  work( 1 ) = 2
218  IF( wantz )
219  $ a( 1, 1 ) = one
220  RETURN
221  END IF
222 *
223 * Get machine constants.
224 *
225  safmin = dlamch( 'Safe minimum' )
226  eps = dlamch( 'Precision' )
227  smlnum = safmin / eps
228  bignum = one / smlnum
229  rmin = sqrt( smlnum )
230  rmax = sqrt( bignum )
231 *
232 * Scale matrix to allowable range, if necessary.
233 *
234  anrm = dlansy( 'M', uplo, n, a, lda, work )
235  iscale = 0
236  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
237  iscale = 1
238  sigma = rmin / anrm
239  ELSE IF( anrm.GT.rmax ) THEN
240  iscale = 1
241  sigma = rmax / anrm
242  END IF
243  IF( iscale.EQ.1 )
244  $ CALL dlascl( uplo, 0, 0, one, sigma, n, n, a, lda, info )
245 *
246 * Call DSYTRD to reduce symmetric matrix to tridiagonal form.
247 *
248  inde = 1
249  indtau = inde + n
250  indwrk = indtau + n
251  llwork = lwork - indwrk + 1
252  CALL dsytrd( uplo, n, a, lda, w, work( inde ), work( indtau ),
253  $ work( indwrk ), llwork, iinfo )
254 *
255 * For eigenvalues only, call DSTERF. For eigenvectors, first call
256 * DORGTR to generate the orthogonal matrix, then call DSTEQR.
257 *
258  IF( .NOT.wantz ) THEN
259  CALL dsterf( n, w, work( inde ), info )
260  ELSE
261  CALL dorgtr( uplo, n, a, lda, work( indtau ), work( indwrk ),
262  $ llwork, iinfo )
263  CALL dsteqr( jobz, n, w, work( inde ), a, lda, work( indtau ),
264  $ info )
265  END IF
266 *
267 * If matrix was scaled, then rescale eigenvalues appropriately.
268 *
269  IF( iscale.EQ.1 ) THEN
270  IF( info.EQ.0 ) THEN
271  imax = n
272  ELSE
273  imax = info - 1
274  END IF
275  CALL dscal( imax, one / sigma, w, 1 )
276  END IF
277 *
278 * Set WORK(1) to optimal workspace size.
279 *
280  work( 1 ) = lwkopt
281 *
282  RETURN
283 *
284 * End of DSYEV
285 *
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
subroutine dorgtr(UPLO, N, A, LDA, TAU, WORK, LWORK, INFO)
DORGTR
Definition: dorgtr.f:125
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 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 dsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
DSTEQR
Definition: dsteqr.f:133
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
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: