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

◆ sspevd()

subroutine sspevd ( character  jobz,
character  uplo,
integer  n,
real, dimension( * )  ap,
real, dimension( * )  w,
real, dimension( ldz, * )  z,
integer  ldz,
real, dimension( * )  work,
integer  lwork,
integer, dimension( * )  iwork,
integer  liwork,
integer  info 
)

SSPEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices

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

Purpose:
 SSPEVD computes all the eigenvalues and, optionally, eigenvectors
 of a real symmetric matrix A in packed storage. If eigenvectors are
 desired, it uses a divide and conquer algorithm.
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]AP
          AP is REAL array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the symmetric matrix
          A, packed columnwise in a linear array.  The j-th column of A
          is stored in the array AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.

          On exit, AP is overwritten by values generated during the
          reduction to tridiagonal form.  If UPLO = 'U', the diagonal
          and first superdiagonal of the tridiagonal matrix T overwrite
          the corresponding elements of A, and if UPLO = 'L', the
          diagonal and first subdiagonal of T overwrite the
          corresponding elements of A.
[out]W
          W is REAL array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is REAL array, dimension (LDZ, N)
          If JOBZ = 'V', then if INFO = 0, Z contains the orthonormal
          eigenvectors of the matrix A, with the i-th column of Z
          holding the eigenvector associated with W(i).
          If JOBZ = 'N', then Z is not referenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= 1, and if
          JOBZ = 'V', LDZ >= max(1,N).
[out]WORK
          WORK is REAL array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the required 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.
          If JOBZ = 'V' and N > 1, LWORK must be at least
                                                 1 + 6*N + N**2.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the required 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 required LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of the array IWORK.
          If JOBZ  = 'N' or 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 required 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, 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.

Definition at line 170 of file sspevd.f.

172*
173* -- LAPACK driver routine --
174* -- LAPACK is a software package provided by Univ. of Tennessee, --
175* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
176*
177* .. Scalar Arguments ..
178 CHARACTER JOBZ, UPLO
179 INTEGER INFO, LDZ, LIWORK, LWORK, N
180* ..
181* .. Array Arguments ..
182 INTEGER IWORK( * )
183 REAL AP( * ), W( * ), WORK( * ), Z( LDZ, * )
184* ..
185*
186* =====================================================================
187*
188* .. Parameters ..
189 REAL ZERO, ONE
190 parameter( zero = 0.0e+0, one = 1.0e+0 )
191* ..
192* .. Local Scalars ..
193 LOGICAL LQUERY, WANTZ
194 INTEGER IINFO, INDE, INDTAU, INDWRK, ISCALE, LIWMIN,
195 $ LLWORK, LWMIN
196 REAL ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
197 $ SMLNUM
198* ..
199* .. External Functions ..
200 LOGICAL LSAME
201 REAL SLAMCH, SLANSP, SROUNDUP_LWORK
203* ..
204* .. External Subroutines ..
205 EXTERNAL sopmtr, sscal, ssptrd, sstedc, ssterf, xerbla
206* ..
207* .. Intrinsic Functions ..
208 INTRINSIC sqrt
209* ..
210* .. Executable Statements ..
211*
212* Test the input parameters.
213*
214 wantz = lsame( jobz, 'V' )
215 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
216*
217 info = 0
218 IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
219 info = -1
220 ELSE IF( .NOT.( lsame( uplo, 'U' ) .OR. lsame( uplo, 'L' ) ) )
221 $ THEN
222 info = -2
223 ELSE IF( n.LT.0 ) THEN
224 info = -3
225 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
226 info = -7
227 END IF
228*
229 IF( info.EQ.0 ) THEN
230 IF( n.LE.1 ) THEN
231 liwmin = 1
232 lwmin = 1
233 ELSE
234 IF( wantz ) THEN
235 liwmin = 3 + 5*n
236 lwmin = 1 + 6*n + n**2
237 ELSE
238 liwmin = 1
239 lwmin = 2*n
240 END IF
241 END IF
242 iwork( 1 ) = liwmin
243 work( 1 ) = sroundup_lwork(lwmin)
244*
245 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
246 info = -9
247 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
248 info = -11
249 END IF
250 END IF
251*
252 IF( info.NE.0 ) THEN
253 CALL xerbla( 'SSPEVD', -info )
254 RETURN
255 ELSE IF( lquery ) THEN
256 RETURN
257 END IF
258*
259* Quick return if possible
260*
261 IF( n.EQ.0 )
262 $ RETURN
263*
264 IF( n.EQ.1 ) THEN
265 w( 1 ) = ap( 1 )
266 IF( wantz )
267 $ z( 1, 1 ) = one
268 RETURN
269 END IF
270*
271* Get machine constants.
272*
273 safmin = slamch( 'Safe minimum' )
274 eps = slamch( 'Precision' )
275 smlnum = safmin / eps
276 bignum = one / smlnum
277 rmin = sqrt( smlnum )
278 rmax = sqrt( bignum )
279*
280* Scale matrix to allowable range, if necessary.
281*
282 anrm = slansp( 'M', uplo, n, ap, work )
283 iscale = 0
284 IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
285 iscale = 1
286 sigma = rmin / anrm
287 ELSE IF( anrm.GT.rmax ) THEN
288 iscale = 1
289 sigma = rmax / anrm
290 END IF
291 IF( iscale.EQ.1 ) THEN
292 CALL sscal( ( n*( n+1 ) ) / 2, sigma, ap, 1 )
293 END IF
294*
295* Call SSPTRD to reduce symmetric packed matrix to tridiagonal form.
296*
297 inde = 1
298 indtau = inde + n
299 CALL ssptrd( uplo, n, ap, w, work( inde ), work( indtau ), iinfo )
300*
301* For eigenvalues only, call SSTERF. For eigenvectors, first call
302* SSTEDC to generate the eigenvector matrix, WORK(INDWRK), of the
303* tridiagonal matrix, then call SOPMTR to multiply it by the
304* Householder transformations represented in AP.
305*
306 IF( .NOT.wantz ) THEN
307 CALL ssterf( n, w, work( inde ), info )
308 ELSE
309 indwrk = indtau + n
310 llwork = lwork - indwrk + 1
311 CALL sstedc( 'I', n, w, work( inde ), z, ldz, work( indwrk ),
312 $ llwork, iwork, liwork, info )
313 CALL sopmtr( 'L', uplo, 'N', n, n, ap, work( indtau ), z, ldz,
314 $ work( indwrk ), iinfo )
315 END IF
316*
317* If matrix was scaled, then rescale eigenvalues appropriately.
318*
319 IF( iscale.EQ.1 )
320 $ CALL sscal( n, one / sigma, w, 1 )
321*
322 work( 1 ) = sroundup_lwork(lwmin)
323 iwork( 1 ) = liwmin
324 RETURN
325*
326* End of SSPEVD
327*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ssptrd(uplo, n, ap, d, e, tau, info)
SSPTRD
Definition ssptrd.f:150
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function slansp(norm, uplo, n, ap, work)
SLANSP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition slansp.f:114
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 sopmtr(side, uplo, trans, m, n, ap, tau, c, ldc, work, info)
SOPMTR
Definition sopmtr.f:150
Here is the call graph for this function:
Here is the caller graph for this function: