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