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

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

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

Purpose:
 SSYEVD 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, SSYEVD needs N**2 more
 workspace than SSYEVX.
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 REAL 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 REAL array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]WORK
          WORK is REAL 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
November 2011
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 185 of file ssyevd.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: