LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
dspev.f
Go to the documentation of this file.
1 *> \brief <b> DSPEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrices</b>
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download DSPEV + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dspev.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dspev.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dspev.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE DSPEV( JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, INFO )
22 *
23 * .. Scalar Arguments ..
24 * CHARACTER JOBZ, UPLO
25 * INTEGER INFO, LDZ, N
26 * ..
27 * .. Array Arguments ..
28 * DOUBLE PRECISION AP( * ), W( * ), WORK( * ), Z( LDZ, * )
29 * ..
30 *
31 *
32 *> \par Purpose:
33 * =============
34 *>
35 *> \verbatim
36 *>
37 *> DSPEV computes all the eigenvalues and, optionally, eigenvectors of a
38 *> real symmetric matrix A in packed storage.
39 *> \endverbatim
40 *
41 * Arguments:
42 * ==========
43 *
44 *> \param[in] JOBZ
45 *> \verbatim
46 *> JOBZ is CHARACTER*1
47 *> = 'N': Compute eigenvalues only;
48 *> = 'V': Compute eigenvalues and eigenvectors.
49 *> \endverbatim
50 *>
51 *> \param[in] UPLO
52 *> \verbatim
53 *> UPLO is CHARACTER*1
54 *> = 'U': Upper triangle of A is stored;
55 *> = 'L': Lower triangle of A is stored.
56 *> \endverbatim
57 *>
58 *> \param[in] N
59 *> \verbatim
60 *> N is INTEGER
61 *> The order of the matrix A. N >= 0.
62 *> \endverbatim
63 *>
64 *> \param[in,out] AP
65 *> \verbatim
66 *> AP is DOUBLE PRECISION array, dimension (N*(N+1)/2)
67 *> On entry, the upper or lower triangle of the symmetric matrix
68 *> A, packed columnwise in a linear array. The j-th column of A
69 *> is stored in the array AP as follows:
70 *> if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
71 *> if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.
72 *>
73 *> On exit, AP is overwritten by values generated during the
74 *> reduction to tridiagonal form. If UPLO = 'U', the diagonal
75 *> and first superdiagonal of the tridiagonal matrix T overwrite
76 *> the corresponding elements of A, and if UPLO = 'L', the
77 *> diagonal and first subdiagonal of T overwrite the
78 *> corresponding elements of A.
79 *> \endverbatim
80 *>
81 *> \param[out] W
82 *> \verbatim
83 *> W is DOUBLE PRECISION array, dimension (N)
84 *> If INFO = 0, the eigenvalues in ascending order.
85 *> \endverbatim
86 *>
87 *> \param[out] Z
88 *> \verbatim
89 *> Z is DOUBLE PRECISION array, dimension (LDZ, N)
90 *> If JOBZ = 'V', then if INFO = 0, Z contains the orthonormal
91 *> eigenvectors of the matrix A, with the i-th column of Z
92 *> holding the eigenvector associated with W(i).
93 *> If JOBZ = 'N', then Z is not referenced.
94 *> \endverbatim
95 *>
96 *> \param[in] LDZ
97 *> \verbatim
98 *> LDZ is INTEGER
99 *> The leading dimension of the array Z. LDZ >= 1, and if
100 *> JOBZ = 'V', LDZ >= max(1,N).
101 *> \endverbatim
102 *>
103 *> \param[out] WORK
104 *> \verbatim
105 *> WORK is DOUBLE PRECISION array, dimension (3*N)
106 *> \endverbatim
107 *>
108 *> \param[out] INFO
109 *> \verbatim
110 *> INFO is INTEGER
111 *> = 0: successful exit.
112 *> < 0: if INFO = -i, the i-th argument had an illegal value.
113 *> > 0: if INFO = i, the algorithm failed to converge; i
114 *> off-diagonal elements of an intermediate tridiagonal
115 *> form did not converge to zero.
116 *> \endverbatim
117 *
118 * Authors:
119 * ========
120 *
121 *> \author Univ. of Tennessee
122 *> \author Univ. of California Berkeley
123 *> \author Univ. of Colorado Denver
124 *> \author NAG Ltd.
125 *
126 *> \date November 2011
127 *
128 *> \ingroup doubleOTHEReigen
129 *
130 * =====================================================================
131  SUBROUTINE dspev( JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, INFO )
132 *
133 * -- LAPACK driver routine (version 3.4.0) --
134 * -- LAPACK is a software package provided by Univ. of Tennessee, --
135 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
136 * November 2011
137 *
138 * .. Scalar Arguments ..
139  CHARACTER JOBZ, UPLO
140  INTEGER INFO, LDZ, N
141 * ..
142 * .. Array Arguments ..
143  DOUBLE PRECISION AP( * ), W( * ), WORK( * ), Z( ldz, * )
144 * ..
145 *
146 * =====================================================================
147 *
148 * .. Parameters ..
149  DOUBLE PRECISION ZERO, ONE
150  parameter ( zero = 0.0d0, one = 1.0d0 )
151 * ..
152 * .. Local Scalars ..
153  LOGICAL WANTZ
154  INTEGER IINFO, IMAX, INDE, INDTAU, INDWRK, ISCALE
155  DOUBLE PRECISION ANRM, BIGNUM, EPS, RMAX, RMIN, SAFMIN, SIGMA,
156  $ smlnum
157 * ..
158 * .. External Functions ..
159  LOGICAL LSAME
160  DOUBLE PRECISION DLAMCH, DLANSP
161  EXTERNAL lsame, dlamch, dlansp
162 * ..
163 * .. External Subroutines ..
164  EXTERNAL dopgtr, dscal, dsptrd, dsteqr, dsterf, xerbla
165 * ..
166 * .. Intrinsic Functions ..
167  INTRINSIC sqrt
168 * ..
169 * .. Executable Statements ..
170 *
171 * Test the input parameters.
172 *
173  wantz = lsame( jobz, 'V' )
174 *
175  info = 0
176  IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
177  info = -1
178  ELSE IF( .NOT.( lsame( uplo, 'U' ) .OR. lsame( uplo, 'L' ) ) )
179  $ THEN
180  info = -2
181  ELSE IF( n.LT.0 ) THEN
182  info = -3
183  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
184  info = -7
185  END IF
186 *
187  IF( info.NE.0 ) THEN
188  CALL xerbla( 'DSPEV ', -info )
189  RETURN
190  END IF
191 *
192 * Quick return if possible
193 *
194  IF( n.EQ.0 )
195  $ RETURN
196 *
197  IF( n.EQ.1 ) THEN
198  w( 1 ) = ap( 1 )
199  IF( wantz )
200  $ z( 1, 1 ) = one
201  RETURN
202  END IF
203 *
204 * Get machine constants.
205 *
206  safmin = dlamch( 'Safe minimum' )
207  eps = dlamch( 'Precision' )
208  smlnum = safmin / eps
209  bignum = one / smlnum
210  rmin = sqrt( smlnum )
211  rmax = sqrt( bignum )
212 *
213 * Scale matrix to allowable range, if necessary.
214 *
215  anrm = dlansp( 'M', uplo, n, ap, work )
216  iscale = 0
217  IF( anrm.GT.zero .AND. anrm.LT.rmin ) THEN
218  iscale = 1
219  sigma = rmin / anrm
220  ELSE IF( anrm.GT.rmax ) THEN
221  iscale = 1
222  sigma = rmax / anrm
223  END IF
224  IF( iscale.EQ.1 ) THEN
225  CALL dscal( ( n*( n+1 ) ) / 2, sigma, ap, 1 )
226  END IF
227 *
228 * Call DSPTRD to reduce symmetric packed matrix to tridiagonal form.
229 *
230  inde = 1
231  indtau = inde + n
232  CALL dsptrd( uplo, n, ap, w, work( inde ), work( indtau ), iinfo )
233 *
234 * For eigenvalues only, call DSTERF. For eigenvectors, first call
235 * DOPGTR to generate the orthogonal matrix, then call DSTEQR.
236 *
237  IF( .NOT.wantz ) THEN
238  CALL dsterf( n, w, work( inde ), info )
239  ELSE
240  indwrk = indtau + n
241  CALL dopgtr( uplo, n, ap, work( indtau ), z, ldz,
242  $ work( indwrk ), iinfo )
243  CALL dsteqr( jobz, n, w, work( inde ), z, ldz, work( indtau ),
244  $ info )
245  END IF
246 *
247 * If matrix was scaled, then rescale eigenvalues appropriately.
248 *
249  IF( iscale.EQ.1 ) THEN
250  IF( info.EQ.0 ) THEN
251  imax = n
252  ELSE
253  imax = info - 1
254  END IF
255  CALL dscal( imax, one / sigma, w, 1 )
256  END IF
257 *
258  RETURN
259 *
260 * End of DSPEV
261 *
262  END
subroutine dsterf(N, D, E, INFO)
DSTERF
Definition: dsterf.f:88
subroutine dsteqr(COMPZ, N, D, E, Z, LDZ, WORK, INFO)
DSTEQR
Definition: dsteqr.f:133
subroutine dspev(JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, INFO)
DSPEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: dspev.f:132
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 dopgtr(UPLO, N, AP, TAU, Q, LDQ, WORK, INFO)
DOPGTR
Definition: dopgtr.f:116