LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dspevd ( character  JOBZ,
character  UPLO,
integer  N,
double precision, dimension( * )  AP,
double precision, dimension( * )  W,
double precision, dimension( ldz, * )  Z,
integer  LDZ,
double precision, dimension( * )  WORK,
integer  LWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

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

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

Purpose:
 DSPEVD 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.

 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.
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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is DOUBLE PRECISION 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 DOUBLE PRECISION array,
                                         dimension (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.
Date
November 2011

Definition at line 181 of file dspevd.f.

181 *
182 * -- LAPACK driver routine (version 3.4.0) --
183 * -- LAPACK is a software package provided by Univ. of Tennessee, --
184 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
185 * November 2011
186 *
187 * .. Scalar Arguments ..
188  CHARACTER jobz, uplo
189  INTEGER info, ldz, liwork, lwork, n
190 * ..
191 * .. Array Arguments ..
192  INTEGER iwork( * )
193  DOUBLE PRECISION ap( * ), w( * ), work( * ), z( ldz, * )
194 * ..
195 *
196 * =====================================================================
197 *
198 * .. Parameters ..
199  DOUBLE PRECISION zero, one
200  parameter ( zero = 0.0d+0, one = 1.0d+0 )
201 * ..
202 * .. Local Scalars ..
203  LOGICAL lquery, wantz
204  INTEGER iinfo, inde, indtau, indwrk, iscale, liwmin,
205  $ llwork, lwmin
206  DOUBLE PRECISION anrm, bignum, eps, rmax, rmin, safmin, sigma,
207  $ smlnum
208 * ..
209 * .. External Functions ..
210  LOGICAL lsame
211  DOUBLE PRECISION dlamch, dlansp
212  EXTERNAL lsame, dlamch, dlansp
213 * ..
214 * .. External Subroutines ..
215  EXTERNAL dopmtr, dscal, dsptrd, dstedc, dsterf, xerbla
216 * ..
217 * .. Intrinsic Functions ..
218  INTRINSIC sqrt
219 * ..
220 * .. Executable Statements ..
221 *
222 * Test the input parameters.
223 *
224  wantz = lsame( jobz, 'V' )
225  lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
226 *
227  info = 0
228  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
229  info = -1
230  ELSE IF( .NOT.( lsame( uplo, 'U' ) .OR. lsame( uplo, 'L' ) ) )
231  $ THEN
232  info = -2
233  ELSE IF( n.LT.0 ) THEN
234  info = -3
235  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
236  info = -7
237  END IF
238 *
239  IF( info.EQ.0 ) THEN
240  IF( n.LE.1 ) THEN
241  liwmin = 1
242  lwmin = 1
243  ELSE
244  IF( wantz ) THEN
245  liwmin = 3 + 5*n
246  lwmin = 1 + 6*n + n**2
247  ELSE
248  liwmin = 1
249  lwmin = 2*n
250  END IF
251  END IF
252  iwork( 1 ) = liwmin
253  work( 1 ) = lwmin
254 *
255  IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
256  info = -9
257  ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
258  info = -11
259  END IF
260  END IF
261 *
262  IF( info.NE.0 ) THEN
263  CALL xerbla( 'DSPEVD', -info )
264  RETURN
265  ELSE IF( lquery ) THEN
266  RETURN
267  END IF
268 *
269 * Quick return if possible
270 *
271  IF( n.EQ.0 )
272  $ RETURN
273 *
274  IF( n.EQ.1 ) THEN
275  w( 1 ) = ap( 1 )
276  IF( wantz )
277  $ z( 1, 1 ) = one
278  RETURN
279  END IF
280 *
281 * Get machine constants.
282 *
283  safmin = dlamch( 'Safe minimum' )
284  eps = dlamch( 'Precision' )
285  smlnum = safmin / eps
286  bignum = one / smlnum
287  rmin = sqrt( smlnum )
288  rmax = sqrt( bignum )
289 *
290 * Scale matrix to allowable range, if necessary.
291 *
292  anrm = dlansp( 'M', uplo, n, ap, work )
293  iscale = 0
294  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
295  iscale = 1
296  sigma = rmin / anrm
297  ELSE IF( anrm.GT.rmax ) THEN
298  iscale = 1
299  sigma = rmax / anrm
300  END IF
301  IF( iscale.EQ.1 ) THEN
302  CALL dscal( ( n*( n+1 ) ) / 2, sigma, ap, 1 )
303  END IF
304 *
305 * Call DSPTRD to reduce symmetric packed matrix to tridiagonal form.
306 *
307  inde = 1
308  indtau = inde + n
309  CALL dsptrd( uplo, n, ap, w, work( inde ), work( indtau ), iinfo )
310 *
311 * For eigenvalues only, call DSTERF. For eigenvectors, first call
312 * DSTEDC to generate the eigenvector matrix, WORK(INDWRK), of the
313 * tridiagonal matrix, then call DOPMTR to multiply it by the
314 * Householder transformations represented in AP.
315 *
316  IF( .NOT.wantz ) THEN
317  CALL dsterf( n, w, work( inde ), info )
318  ELSE
319  indwrk = indtau + n
320  llwork = lwork - indwrk + 1
321  CALL dstedc( 'I', n, w, work( inde ), z, ldz, work( indwrk ),
322  $ llwork, iwork, liwork, info )
323  CALL dopmtr( 'L', uplo, 'N', n, n, ap, work( indtau ), z, ldz,
324  $ work( indwrk ), iinfo )
325  END IF
326 *
327 * If matrix was scaled, then rescale eigenvalues appropriately.
328 *
329  IF( iscale.EQ.1 )
330  $ CALL dscal( n, one / sigma, w, 1 )
331 *
332  work( 1 ) = lwmin
333  iwork( 1 ) = liwmin
334  RETURN
335 *
336 * End of DSPEVD
337 *
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dopmtr(SIDE, UPLO, TRANS, M, N, AP, TAU, C, LDC, WORK, INFO)
DOPMTR
Definition: dopmtr.f:152
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
subroutine dsptrd(UPLO, N, AP, D, E, TAU, INFO)
DSPTRD
Definition: dsptrd.f:152
subroutine dstedc(COMPZ, N, D, E, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
DSTEDC
Definition: dstedc.f:191
double precision function dlansp(NORM, UPLO, N, AP, WORK)
DLANSP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a symmetric matrix supplied in packed form.
Definition: dlansp.f:116
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: