LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dpteqr()

subroutine dpteqr ( character  compz,
integer  n,
double precision, dimension( * )  d,
double precision, dimension( * )  e,
double precision, dimension( ldz, * )  z,
integer  ldz,
double precision, dimension( * )  work,
integer  info 
)

DPTEQR

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

Purpose:
 DPTEQR computes all eigenvalues and, optionally, eigenvectors of a
 symmetric positive definite tridiagonal matrix by first factoring the
 matrix using DPTTRF, and then calling DBDSQR to compute the singular
 values of the bidiagonal factor.

 This routine computes the eigenvalues of the positive definite
 tridiagonal matrix to high relative accuracy.  This means that if the
 eigenvalues range over many orders of magnitude in size, then the
 small eigenvalues and corresponding eigenvectors will be computed
 more accurately than, for example, with the standard QR method.

 The eigenvectors of a full or band symmetric positive definite matrix
 can also be found if DSYTRD, DSPTRD, or DSBTRD has been used to
 reduce this matrix to tridiagonal form. (The reduction to tridiagonal
 form, however, may preclude the possibility of obtaining high
 relative accuracy in the small eigenvalues of the original matrix, if
 these eigenvalues range over many orders of magnitude.)
Parameters
[in]COMPZ
          COMPZ is CHARACTER*1
          = 'N':  Compute eigenvalues only.
          = 'V':  Compute eigenvectors of original symmetric
                  matrix also.  Array Z contains the orthogonal
                  matrix used to reduce the original matrix to
                  tridiagonal form.
          = 'I':  Compute eigenvectors of tridiagonal matrix also.
[in]N
          N is INTEGER
          The order of the matrix.  N >= 0.
[in,out]D
          D is DOUBLE PRECISION array, dimension (N)
          On entry, the n diagonal elements of the tridiagonal
          matrix.
          On normal exit, D contains the eigenvalues, in descending
          order.
[in,out]E
          E is DOUBLE PRECISION array, dimension (N-1)
          On entry, the (n-1) subdiagonal elements of the tridiagonal
          matrix.
          On exit, E has been destroyed.
[in,out]Z
          Z is DOUBLE PRECISION array, dimension (LDZ, N)
          On entry, if COMPZ = 'V', the orthogonal matrix used in the
          reduction to tridiagonal form.
          On exit, if COMPZ = 'V', the orthonormal eigenvectors of the
          original symmetric matrix;
          if COMPZ = 'I', the orthonormal eigenvectors of the
          tridiagonal matrix.
          If INFO > 0 on exit, Z contains the eigenvectors associated
          with only the stored eigenvalues.
          If  COMPZ = 'N', then Z is not referenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= 1, and if
          COMPZ = 'V' or 'I', LDZ >= max(1,N).
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (4*N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.
          > 0:  if INFO = i, and i is:
                <= N  the Cholesky factorization of the matrix could
                      not be performed because the leading principal
                      minor of order i was not positive.
                > N   the SVD algorithm failed to converge;
                      if INFO = N+i, i off-diagonal elements of the
                      bidiagonal factor did not converge to zero.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 144 of file dpteqr.f.

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*
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:110
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dpttrf(n, d, e, info)
DPTTRF
Definition dpttrf.f:91
Here is the call graph for this function:
Here is the caller graph for this function: