LAPACK 3.3.0

dpteqr.f

Go to the documentation of this file.
00001       SUBROUTINE DPTEQR( COMPZ, N, D, E, Z, LDZ, WORK, INFO )
00002 *
00003 *  -- LAPACK routine (version 3.2) --
00004 *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
00005 *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
00006 *     November 2006
00007 *
00008 *     .. Scalar Arguments ..
00009       CHARACTER          COMPZ
00010       INTEGER            INFO, LDZ, N
00011 *     ..
00012 *     .. Array Arguments ..
00013       DOUBLE PRECISION   D( * ), E( * ), WORK( * ), Z( LDZ, * )
00014 *     ..
00015 *
00016 *  Purpose
00017 *  =======
00018 *
00019 *  DPTEQR computes all eigenvalues and, optionally, eigenvectors of a
00020 *  symmetric positive definite tridiagonal matrix by first factoring the
00021 *  matrix using DPTTRF, and then calling DBDSQR to compute the singular
00022 *  values of the bidiagonal factor.
00023 *
00024 *  This routine computes the eigenvalues of the positive definite
00025 *  tridiagonal matrix to high relative accuracy.  This means that if the
00026 *  eigenvalues range over many orders of magnitude in size, then the
00027 *  small eigenvalues and corresponding eigenvectors will be computed
00028 *  more accurately than, for example, with the standard QR method.
00029 *
00030 *  The eigenvectors of a full or band symmetric positive definite matrix
00031 *  can also be found if DSYTRD, DSPTRD, or DSBTRD has been used to
00032 *  reduce this matrix to tridiagonal form. (The reduction to tridiagonal
00033 *  form, however, may preclude the possibility of obtaining high
00034 *  relative accuracy in the small eigenvalues of the original matrix, if
00035 *  these eigenvalues range over many orders of magnitude.)
00036 *
00037 *  Arguments
00038 *  =========
00039 *
00040 *  COMPZ   (input) CHARACTER*1
00041 *          = 'N':  Compute eigenvalues only.
00042 *          = 'V':  Compute eigenvectors of original symmetric
00043 *                  matrix also.  Array Z contains the orthogonal
00044 *                  matrix used to reduce the original matrix to
00045 *                  tridiagonal form.
00046 *          = 'I':  Compute eigenvectors of tridiagonal matrix also.
00047 *
00048 *  N       (input) INTEGER
00049 *          The order of the matrix.  N >= 0.
00050 *
00051 *  D       (input/output) DOUBLE PRECISION array, dimension (N)
00052 *          On entry, the n diagonal elements of the tridiagonal
00053 *          matrix.
00054 *          On normal exit, D contains the eigenvalues, in descending
00055 *          order.
00056 *
00057 *  E       (input/output) DOUBLE PRECISION array, dimension (N-1)
00058 *          On entry, the (n-1) subdiagonal elements of the tridiagonal
00059 *          matrix.
00060 *          On exit, E has been destroyed.
00061 *
00062 *  Z       (input/output) DOUBLE PRECISION array, dimension (LDZ, N)
00063 *          On entry, if COMPZ = 'V', the orthogonal matrix used in the
00064 *          reduction to tridiagonal form.
00065 *          On exit, if COMPZ = 'V', the orthonormal eigenvectors of the
00066 *          original symmetric matrix;
00067 *          if COMPZ = 'I', the orthonormal eigenvectors of the
00068 *          tridiagonal matrix.
00069 *          If INFO > 0 on exit, Z contains the eigenvectors associated
00070 *          with only the stored eigenvalues.
00071 *          If  COMPZ = 'N', then Z is not referenced.
00072 *
00073 *  LDZ     (input) INTEGER
00074 *          The leading dimension of the array Z.  LDZ >= 1, and if
00075 *          COMPZ = 'V' or 'I', LDZ >= max(1,N).
00076 *
00077 *  WORK    (workspace) DOUBLE PRECISION array, dimension (4*N)
00078 *
00079 *  INFO    (output) INTEGER
00080 *          = 0:  successful exit.
00081 *          < 0:  if INFO = -i, the i-th argument had an illegal value.
00082 *          > 0:  if INFO = i, and i is:
00083 *                <= N  the Cholesky factorization of the matrix could
00084 *                      not be performed because the i-th principal minor
00085 *                      was not positive definite.
00086 *                > N   the SVD algorithm failed to converge;
00087 *                      if INFO = N+i, i off-diagonal elements of the
00088 *                      bidiagonal factor did not converge to zero.
00089 *
00090 *  =====================================================================
00091 *
00092 *     .. Parameters ..
00093       DOUBLE PRECISION   ZERO, ONE
00094       PARAMETER          ( ZERO = 0.0D0, ONE = 1.0D0 )
00095 *     ..
00096 *     .. External Functions ..
00097       LOGICAL            LSAME
00098       EXTERNAL           LSAME
00099 *     ..
00100 *     .. External Subroutines ..
00101       EXTERNAL           DBDSQR, DLASET, DPTTRF, XERBLA
00102 *     ..
00103 *     .. Local Arrays ..
00104       DOUBLE PRECISION   C( 1, 1 ), VT( 1, 1 )
00105 *     ..
00106 *     .. Local Scalars ..
00107       INTEGER            I, ICOMPZ, NRU
00108 *     ..
00109 *     .. Intrinsic Functions ..
00110       INTRINSIC          MAX, SQRT
00111 *     ..
00112 *     .. Executable Statements ..
00113 *
00114 *     Test the input parameters.
00115 *
00116       INFO = 0
00117 *
00118       IF( LSAME( COMPZ, 'N' ) ) THEN
00119          ICOMPZ = 0
00120       ELSE IF( LSAME( COMPZ, 'V' ) ) THEN
00121          ICOMPZ = 1
00122       ELSE IF( LSAME( COMPZ, 'I' ) ) THEN
00123          ICOMPZ = 2
00124       ELSE
00125          ICOMPZ = -1
00126       END IF
00127       IF( ICOMPZ.LT.0 ) THEN
00128          INFO = -1
00129       ELSE IF( N.LT.0 ) THEN
00130          INFO = -2
00131       ELSE IF( ( LDZ.LT.1 ) .OR. ( ICOMPZ.GT.0 .AND. LDZ.LT.MAX( 1,
00132      $         N ) ) ) THEN
00133          INFO = -6
00134       END IF
00135       IF( INFO.NE.0 ) THEN
00136          CALL XERBLA( 'DPTEQR', -INFO )
00137          RETURN
00138       END IF
00139 *
00140 *     Quick return if possible
00141 *
00142       IF( N.EQ.0 )
00143      $   RETURN
00144 *
00145       IF( N.EQ.1 ) THEN
00146          IF( ICOMPZ.GT.0 )
00147      $      Z( 1, 1 ) = ONE
00148          RETURN
00149       END IF
00150       IF( ICOMPZ.EQ.2 )
00151      $   CALL DLASET( 'Full', N, N, ZERO, ONE, Z, LDZ )
00152 *
00153 *     Call DPTTRF to factor the matrix.
00154 *
00155       CALL DPTTRF( N, D, E, INFO )
00156       IF( INFO.NE.0 )
00157      $   RETURN
00158       DO 10 I = 1, N
00159          D( I ) = SQRT( D( I ) )
00160    10 CONTINUE
00161       DO 20 I = 1, N - 1
00162          E( I ) = E( I )*D( I )
00163    20 CONTINUE
00164 *
00165 *     Call DBDSQR to compute the singular values/vectors of the
00166 *     bidiagonal factor.
00167 *
00168       IF( ICOMPZ.GT.0 ) THEN
00169          NRU = N
00170       ELSE
00171          NRU = 0
00172       END IF
00173       CALL DBDSQR( 'Lower', N, 0, NRU, 0, D, E, VT, 1, Z, LDZ, C, 1,
00174      $             WORK, INFO )
00175 *
00176 *     Square the singular values.
00177 *
00178       IF( INFO.EQ.0 ) THEN
00179          DO 30 I = 1, N
00180             D( I ) = D( I )*D( I )
00181    30    CONTINUE
00182       ELSE
00183          INFO = N + INFO
00184       END IF
00185 *
00186       RETURN
00187 *
00188 *     End of DPTEQR
00189 *
00190       END
 All Files Functions