LAPACK  3.9.0
LAPACK: Linear Algebra PACKage

◆ dspevd()

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 (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.
Date
June 2017

Definition at line 180 of file dspevd.f.

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