LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dpteqr.f
Go to the documentation of this file.
1*> \brief \b DPTEQR
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DPTEQR + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dpteqr.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dpteqr.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dpteqr.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DPTEQR( COMPZ, N, D, E, Z, LDZ, WORK, INFO )
20*
21* .. Scalar Arguments ..
22* CHARACTER COMPZ
23* INTEGER INFO, LDZ, N
24* ..
25* .. Array Arguments ..
26* DOUBLE PRECISION D( * ), E( * ), WORK( * ), Z( LDZ, * )
27* ..
28*
29*
30*> \par Purpose:
31* =============
32*>
33*> \verbatim
34*>
35*> DPTEQR computes all eigenvalues and, optionally, eigenvectors of a
36*> symmetric positive definite tridiagonal matrix by first factoring the
37*> matrix using DPTTRF, and then calling DBDSQR to compute the singular
38*> values of the bidiagonal factor.
39*>
40*> This routine computes the eigenvalues of the positive definite
41*> tridiagonal matrix to high relative accuracy. This means that if the
42*> eigenvalues range over many orders of magnitude in size, then the
43*> small eigenvalues and corresponding eigenvectors will be computed
44*> more accurately than, for example, with the standard QR method.
45*>
46*> The eigenvectors of a full or band symmetric positive definite matrix
47*> can also be found if DSYTRD, DSPTRD, or DSBTRD has been used to
48*> reduce this matrix to tridiagonal form. (The reduction to tridiagonal
49*> form, however, may preclude the possibility of obtaining high
50*> relative accuracy in the small eigenvalues of the original matrix, if
51*> these eigenvalues range over many orders of magnitude.)
52*> \endverbatim
53*
54* Arguments:
55* ==========
56*
57*> \param[in] COMPZ
58*> \verbatim
59*> COMPZ is CHARACTER*1
60*> = 'N': Compute eigenvalues only.
61*> = 'V': Compute eigenvectors of original symmetric
62*> matrix also. Array Z contains the orthogonal
63*> matrix used to reduce the original matrix to
64*> tridiagonal form.
65*> = 'I': Compute eigenvectors of tridiagonal matrix also.
66*> \endverbatim
67*>
68*> \param[in] N
69*> \verbatim
70*> N is INTEGER
71*> The order of the matrix. N >= 0.
72*> \endverbatim
73*>
74*> \param[in,out] D
75*> \verbatim
76*> D is DOUBLE PRECISION array, dimension (N)
77*> On entry, the n diagonal elements of the tridiagonal
78*> matrix.
79*> On normal exit, D contains the eigenvalues, in descending
80*> order.
81*> \endverbatim
82*>
83*> \param[in,out] E
84*> \verbatim
85*> E is DOUBLE PRECISION array, dimension (N-1)
86*> On entry, the (n-1) subdiagonal elements of the tridiagonal
87*> matrix.
88*> On exit, E has been destroyed.
89*> \endverbatim
90*>
91*> \param[in,out] Z
92*> \verbatim
93*> Z is DOUBLE PRECISION array, dimension (LDZ, N)
94*> On entry, if COMPZ = 'V', the orthogonal matrix used in the
95*> reduction to tridiagonal form.
96*> On exit, if COMPZ = 'V', the orthonormal eigenvectors of the
97*> original symmetric matrix;
98*> if COMPZ = 'I', the orthonormal eigenvectors of the
99*> tridiagonal matrix.
100*> If INFO > 0 on exit, Z contains the eigenvectors associated
101*> with only the stored eigenvalues.
102*> If COMPZ = 'N', then Z is not referenced.
103*> \endverbatim
104*>
105*> \param[in] LDZ
106*> \verbatim
107*> LDZ is INTEGER
108*> The leading dimension of the array Z. LDZ >= 1, and if
109*> COMPZ = 'V' or 'I', LDZ >= max(1,N).
110*> \endverbatim
111*>
112*> \param[out] WORK
113*> \verbatim
114*> WORK is DOUBLE PRECISION array, dimension (4*N)
115*> \endverbatim
116*>
117*> \param[out] INFO
118*> \verbatim
119*> INFO is INTEGER
120*> = 0: successful exit.
121*> < 0: if INFO = -i, the i-th argument had an illegal value.
122*> > 0: if INFO = i, and i is:
123*> <= N the Cholesky factorization of the matrix could
124*> not be performed because the leading principal
125*> minor of order i was not positive.
126*> > N the SVD algorithm failed to converge;
127*> if INFO = N+i, i off-diagonal elements of the
128*> bidiagonal factor did not converge to zero.
129*> \endverbatim
130*
131* Authors:
132* ========
133*
134*> \author Univ. of Tennessee
135*> \author Univ. of California Berkeley
136*> \author Univ. of Colorado Denver
137*> \author NAG Ltd.
138*
139*> \ingroup pteqr
140*
141* =====================================================================
142 SUBROUTINE dpteqr( COMPZ, N, D, E, Z, LDZ, WORK, INFO )
143*
144* -- LAPACK computational routine --
145* -- LAPACK is a software package provided by Univ. of Tennessee, --
146* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
147*
148* .. Scalar Arguments ..
149 CHARACTER COMPZ
150 INTEGER INFO, LDZ, N
151* ..
152* .. Array Arguments ..
153 DOUBLE PRECISION D( * ), E( * ), WORK( * ), Z( LDZ, * )
154* ..
155*
156* =====================================================================
157*
158* .. Parameters ..
159 DOUBLE PRECISION ZERO, ONE
160 parameter( zero = 0.0d0, one = 1.0d0 )
161* ..
162* .. External Functions ..
163 LOGICAL LSAME
164 EXTERNAL lsame
165* ..
166* .. External Subroutines ..
167 EXTERNAL dbdsqr, dlaset, dpttrf, xerbla
168* ..
169* .. Local Arrays ..
170 DOUBLE PRECISION C( 1, 1 ), VT( 1, 1 )
171* ..
172* .. Local Scalars ..
173 INTEGER I, ICOMPZ, NRU
174* ..
175* .. Intrinsic Functions ..
176 INTRINSIC max, sqrt
177* ..
178* .. Executable Statements ..
179*
180* Test the input parameters.
181*
182 info = 0
183*
184 IF( lsame( compz, 'N' ) ) THEN
185 icompz = 0
186 ELSE IF( lsame( compz, 'V' ) ) THEN
187 icompz = 1
188 ELSE IF( lsame( compz, 'I' ) ) THEN
189 icompz = 2
190 ELSE
191 icompz = -1
192 END IF
193 IF( icompz.LT.0 ) THEN
194 info = -1
195 ELSE IF( n.LT.0 ) THEN
196 info = -2
197 ELSE IF( ( ldz.LT.1 ) .OR. ( icompz.GT.0 .AND. ldz.LT.max( 1,
198 $ n ) ) ) THEN
199 info = -6
200 END IF
201 IF( info.NE.0 ) THEN
202 CALL xerbla( 'DPTEQR', -info )
203 RETURN
204 END IF
205*
206* Quick return if possible
207*
208 IF( n.EQ.0 )
209 $ RETURN
210*
211 IF( n.EQ.1 ) THEN
212 IF( icompz.GT.0 )
213 $ z( 1, 1 ) = one
214 RETURN
215 END IF
216 IF( icompz.EQ.2 )
217 $ CALL dlaset( 'Full', n, n, zero, one, z, ldz )
218*
219* Call DPTTRF to factor the matrix.
220*
221 CALL dpttrf( n, d, e, info )
222 IF( info.NE.0 )
223 $ RETURN
224 DO 10 i = 1, n
225 d( i ) = sqrt( d( i ) )
226 10 CONTINUE
227 DO 20 i = 1, n - 1
228 e( i ) = e( i )*d( i )
229 20 CONTINUE
230*
231* Call DBDSQR to compute the singular values/vectors of the
232* bidiagonal factor.
233*
234 IF( icompz.GT.0 ) THEN
235 nru = n
236 ELSE
237 nru = 0
238 END IF
239 CALL dbdsqr( 'Lower', n, 0, nru, 0, d, e, vt, 1, z, ldz, c, 1,
240 $ work, info )
241*
242* Square the singular values.
243*
244 IF( info.EQ.0 ) THEN
245 DO 30 i = 1, n
246 d( i ) = d( i )*d( i )
247 30 CONTINUE
248 ELSE
249 info = n + info
250 END IF
251*
252 RETURN
253*
254* End of DPTEQR
255*
256 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dbdsqr(uplo, n, ncvt, nru, ncc, d, e, vt, ldvt, u, ldu, c, ldc, work, info)
DBDSQR
Definition dbdsqr.f:241
subroutine dlaset(uplo, m, n, alpha, beta, a, lda)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition dlaset.f:108
subroutine dpteqr(compz, n, d, e, z, ldz, work, info)
DPTEQR
Definition dpteqr.f:143
subroutine dpttrf(n, d, e, info)
DPTTRF
Definition dpttrf.f:89