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

◆ ssyevd()

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.

 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.
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 174 of file ssyevd.f.

176*
177* -- LAPACK driver routine --
178* -- LAPACK is a software package provided by Univ. of Tennessee, --
179* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
180*
181* .. Scalar Arguments ..
182 CHARACTER JOBZ, UPLO
183 INTEGER INFO, LDA, LIWORK, LWORK, N
184* ..
185* .. Array Arguments ..
186 INTEGER IWORK( * )
187 REAL A( LDA, * ), W( * ), WORK( * )
188* ..
189*
190* =====================================================================
191*
192* .. Parameters ..
193 REAL ZERO, ONE
194 parameter( zero = 0.0e+0, one = 1.0e+0 )
195* ..
196* .. Local Scalars ..
197*
198 LOGICAL LOWER, LQUERY, WANTZ
199 INTEGER IINFO, INDE, INDTAU, INDWK2, INDWRK, ISCALE,
200 $ LIOPT, LIWMIN, LLWORK, LLWRK2, LOPT, LWMIN
201 REAL ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
202 $ SMLNUM
203* ..
204* .. External Functions ..
205 LOGICAL LSAME
206 INTEGER ILAENV
207 REAL SLAMCH, SLANSY, SROUNDUP_LWORK
209* ..
210* .. External Subroutines ..
211 EXTERNAL slacpy, slascl, sormtr, sscal, sstedc, ssterf,
212 $ ssytrd, xerbla
213* ..
214* .. Intrinsic Functions ..
215 INTRINSIC max, sqrt
216* ..
217* .. Executable Statements ..
218*
219* Test the input parameters.
220*
221 wantz = lsame( jobz, 'V' )
222 lower = lsame( uplo, 'L' )
223 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
224*
225 info = 0
226 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
227 info = -1
228 ELSE IF( .NOT.( lower .OR. lsame( uplo, 'U' ) ) ) THEN
229 info = -2
230 ELSE IF( n.LT.0 ) THEN
231 info = -3
232 ELSE IF( lda.LT.max( 1, n ) ) THEN
233 info = -5
234 END IF
235*
236 IF( info.EQ.0 ) THEN
237 IF( n.LE.1 ) THEN
238 liwmin = 1
239 lwmin = 1
240 lopt = lwmin
241 liopt = liwmin
242 ELSE
243 IF( wantz ) THEN
244 liwmin = 3 + 5*n
245 lwmin = 1 + 6*n + 2*n**2
246 ELSE
247 liwmin = 1
248 lwmin = 2*n + 1
249 END IF
250 lopt = max( lwmin, 2*n +
251 $ n*ilaenv( 1, 'SSYTRD', uplo, n, -1, -1, -1 ) )
252 liopt = liwmin
253 END IF
254 work( 1 ) = sroundup_lwork(lopt)
255 iwork( 1 ) = liopt
256*
257 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
258 info = -8
259 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
260 info = -10
261 END IF
262 END IF
263*
264 IF( info.NE.0 ) THEN
265 CALL xerbla( 'SSYEVD', -info )
266 RETURN
267 ELSE IF( lquery ) THEN
268 RETURN
269 END IF
270*
271* Quick return if possible
272*
273 IF( n.EQ.0 )
274 $ RETURN
275*
276 IF( n.EQ.1 ) THEN
277 w( 1 ) = a( 1, 1 )
278 IF( wantz )
279 $ a( 1, 1 ) = one
280 RETURN
281 END IF
282*
283* Get machine constants.
284*
285 safmin = slamch( 'Safe minimum' )
286 eps = slamch( 'Precision' )
287 smlnum = safmin / eps
288 bignum = one / smlnum
289 rmin = sqrt( smlnum )
290 rmax = sqrt( bignum )
291*
292* Scale matrix to allowable range, if necessary.
293*
294 anrm = slansy( 'M', uplo, n, a, lda, work )
295 iscale = 0
296 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
297 iscale = 1
298 sigma = rmin / anrm
299 ELSE IF( anrm.GT.rmax ) THEN
300 iscale = 1
301 sigma = rmax / anrm
302 END IF
303 IF( iscale.EQ.1 )
304 $ CALL slascl( uplo, 0, 0, one, sigma, n, n, a, lda, info )
305*
306* Call SSYTRD to reduce symmetric matrix to tridiagonal form.
307*
308 inde = 1
309 indtau = inde + n
310 indwrk = indtau + n
311 llwork = lwork - indwrk + 1
312 indwk2 = indwrk + n*n
313 llwrk2 = lwork - indwk2 + 1
314*
315 CALL ssytrd( uplo, n, a, lda, w, work( inde ), work( indtau ),
316 $ work( indwrk ), llwork, iinfo )
317*
318* For eigenvalues only, call SSTERF. For eigenvectors, first call
319* SSTEDC to generate the eigenvector matrix, WORK(INDWRK), of the
320* tridiagonal matrix, then call SORMTR to multiply it by the
321* Householder transformations stored in A.
322*
323 IF( .NOT.wantz ) THEN
324 CALL ssterf( n, w, work( inde ), info )
325 ELSE
326 CALL sstedc( 'I', n, w, work( inde ), work( indwrk ), n,
327 $ work( indwk2 ), llwrk2, iwork, liwork, info )
328 CALL sormtr( 'L', uplo, 'N', n, n, a, lda, work( indtau ),
329 $ work( indwrk ), n, work( indwk2 ), llwrk2, iinfo )
330 CALL slacpy( 'A', n, n, work( indwrk ), n, a, lda )
331 END IF
332*
333* If matrix was scaled, then rescale eigenvalues appropriately.
334*
335 IF( iscale.EQ.1 )
336 $ CALL sscal( n, one / sigma, w, 1 )
337*
338 work( 1 ) = sroundup_lwork(lopt)
339 iwork( 1 ) = liopt
340*
341 RETURN
342*
343* End of SSYEVD
344*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ssytrd(uplo, n, a, lda, d, e, tau, work, lwork, info)
SSYTRD
Definition ssytrd.f:192
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
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,...
Definition slansy.f:122
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:143
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine sstedc(compz, n, d, e, z, ldz, work, lwork, iwork, liwork, info)
SSTEDC
Definition sstedc.f:182
subroutine ssterf(n, d, e, info)
SSTERF
Definition ssterf.f:86
subroutine sormtr(side, uplo, trans, m, n, a, lda, tau, c, ldc, work, lwork, info)
SORMTR
Definition sormtr.f:172
Here is the call graph for this function:
Here is the caller graph for this function: