LAPACK 3.12.1
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 142 of file dpteqr.f.

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