LAPACK 3.12.1
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 (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 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 171 of file ssyevd.f.

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