1*> \brief \b DSPTRD
2*
3*  =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6*            http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download DSPTRD + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsptrd.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsptrd.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsptrd.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18*  Definition:
19*  ===========
20*
21*       SUBROUTINE DSPTRD( UPLO, N, AP, D, E, TAU, INFO )
22*
23*       .. Scalar Arguments ..
24*       CHARACTER          UPLO
25*       INTEGER            INFO, N
26*       ..
27*       .. Array Arguments ..
28*       DOUBLE PRECISION   AP( * ), D( * ), E( * ), TAU( * )
29*       ..
30*
31*
32*> \par Purpose:
33*  =============
34*>
35*> \verbatim
36*>
37*> DSPTRD reduces a real symmetric matrix A stored in packed form to
38*> symmetric tridiagonal form T by an orthogonal similarity
39*> transformation: Q**T * A * Q = T.
40*> \endverbatim
41*
42*  Arguments:
43*  ==========
44*
45*> \param[in] UPLO
46*> \verbatim
47*>          UPLO is CHARACTER*1
48*>          = 'U':  Upper triangle of A is stored;
49*>          = 'L':  Lower triangle of A is stored.
50*> \endverbatim
51*>
52*> \param[in] N
53*> \verbatim
54*>          N is INTEGER
55*>          The order of the matrix A.  N >= 0.
56*> \endverbatim
57*>
58*> \param[in,out] AP
59*> \verbatim
60*>          AP is DOUBLE PRECISION array, dimension (N*(N+1)/2)
61*>          On entry, the upper or lower triangle of the symmetric matrix
62*>          A, packed columnwise in a linear array.  The j-th column of A
63*>          is stored in the array AP as follows:
64*>          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
65*>          if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.
66*>          On exit, if UPLO = 'U', the diagonal and first superdiagonal
67*>          of A are overwritten by the corresponding elements of the
68*>          tridiagonal matrix T, and the elements above the first
69*>          superdiagonal, with the array TAU, represent the orthogonal
70*>          matrix Q as a product of elementary reflectors; if UPLO
71*>          = 'L', the diagonal and first subdiagonal of A are over-
72*>          written by the corresponding elements of the tridiagonal
73*>          matrix T, and the elements below the first subdiagonal, with
74*>          the array TAU, represent the orthogonal matrix Q as a product
75*>          of elementary reflectors. See Further Details.
76*> \endverbatim
77*>
78*> \param[out] D
79*> \verbatim
80*>          D is DOUBLE PRECISION array, dimension (N)
81*>          The diagonal elements of the tridiagonal matrix T:
82*>          D(i) = A(i,i).
83*> \endverbatim
84*>
85*> \param[out] E
86*> \verbatim
87*>          E is DOUBLE PRECISION array, dimension (N-1)
88*>          The off-diagonal elements of the tridiagonal matrix T:
89*>          E(i) = A(i,i+1) if UPLO = 'U', E(i) = A(i+1,i) if UPLO = 'L'.
90*> \endverbatim
91*>
92*> \param[out] TAU
93*> \verbatim
94*>          TAU is DOUBLE PRECISION array, dimension (N-1)
95*>          The scalar factors of the elementary reflectors (see Further
96*>          Details).
97*> \endverbatim
98*>
99*> \param[out] INFO
100*> \verbatim
101*>          INFO is INTEGER
102*>          = 0:  successful exit
103*>          < 0:  if INFO = -i, the i-th argument had an illegal value
104*> \endverbatim
105*
106*  Authors:
107*  ========
108*
109*> \author Univ. of Tennessee
110*> \author Univ. of California Berkeley
111*> \author Univ. of Colorado Denver
112*> \author NAG Ltd.
113*
114*> \ingroup doubleOTHERcomputational
115*
116*> \par Further Details:
117*  =====================
118*>
119*> \verbatim
120*>
121*>  If UPLO = 'U', the matrix Q is represented as a product of elementary
122*>  reflectors
123*>
124*>     Q = H(n-1) . . . H(2) H(1).
125*>
126*>  Each H(i) has the form
127*>
128*>     H(i) = I - tau * v * v**T
129*>
130*>  where tau is a real scalar, and v is a real vector with
131*>  v(i+1:n) = 0 and v(i) = 1; v(1:i-1) is stored on exit in AP,
132*>  overwriting A(1:i-1,i+1), and tau is stored in TAU(i).
133*>
134*>  If UPLO = 'L', the matrix Q is represented as a product of elementary
135*>  reflectors
136*>
137*>     Q = H(1) H(2) . . . H(n-1).
138*>
139*>  Each H(i) has the form
140*>
141*>     H(i) = I - tau * v * v**T
142*>
143*>  where tau is a real scalar, and v is a real vector with
144*>  v(1:i) = 0 and v(i+1) = 1; v(i+2:n) is stored on exit in AP,
145*>  overwriting A(i+2:n,i), and tau is stored in TAU(i).
146*> \endverbatim
147*>
148*  =====================================================================
149      SUBROUTINE DSPTRD( UPLO, N, AP, D, E, TAU, INFO )
150*
151*  -- LAPACK computational routine --
152*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
153*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
154*
155*     .. Scalar Arguments ..
156      CHARACTER          UPLO
157      INTEGER            INFO, N
158*     ..
159*     .. Array Arguments ..
160      DOUBLE PRECISION   AP( * ), D( * ), E( * ), TAU( * )
161*     ..
162*
163*  =====================================================================
164*
165*     .. Parameters ..
166      DOUBLE PRECISION   ONE, ZERO, HALF
167      PARAMETER          ( ONE = 1.0D0, ZERO = 0.0D0,
168     $                   HALF = 1.0D0 / 2.0D0 )
169*     ..
170*     .. Local Scalars ..
171      LOGICAL            UPPER
172      INTEGER            I, I1, I1I1, II
173      DOUBLE PRECISION   ALPHA, TAUI
174*     ..
175*     .. External Subroutines ..
176      EXTERNAL           DAXPY, DLARFG, DSPMV, DSPR2, XERBLA
177*     ..
178*     .. External Functions ..
179      LOGICAL            LSAME
180      DOUBLE PRECISION   DDOT
181      EXTERNAL           LSAME, DDOT
182*     ..
183*     .. Executable Statements ..
184*
185*     Test the input parameters
186*
187      INFO = 0
188      UPPER = LSAME( UPLO, 'U' )
189      IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
190         INFO = -1
191      ELSE IF( N.LT.0 ) THEN
192         INFO = -2
193      END IF
194      IF( INFO.NE.0 ) THEN
195         CALL XERBLA( 'DSPTRD', -INFO )
196         RETURN
197      END IF
198*
199*     Quick return if possible
200*
201      IF( N.LE.0 )
202     $   RETURN
203*
204      IF( UPPER ) THEN
205*
206*        Reduce the upper triangle of A.
207*        I1 is the index in AP of A(1,I+1).
208*
209         I1 = N*( N-1 ) / 2 + 1
210         DO 10 I = N - 1, 1, -1
211*
212*           Generate elementary reflector H(i) = I - tau * v * v**T
213*           to annihilate A(1:i-1,i+1)
214*
215            CALL DLARFG( I, AP( I1+I-1 ), AP( I1 ), 1, TAUI )
216            E( I ) = AP( I1+I-1 )
217*
218            IF( TAUI.NE.ZERO ) THEN
219*
220*              Apply H(i) from both sides to A(1:i,1:i)
221*
222               AP( I1+I-1 ) = ONE
223*
224*              Compute  y := tau * A * v  storing y in TAU(1:i)
225*
226               CALL DSPMV( UPLO, I, TAUI, AP, AP( I1 ), 1, ZERO, TAU,
227     $                     1 )
228*
229*              Compute  w := y - 1/2 * tau * (y**T *v) * v
230*
231               ALPHA = -HALF*TAUI*DDOT( I, TAU, 1, AP( I1 ), 1 )
232               CALL DAXPY( I, ALPHA, AP( I1 ), 1, TAU, 1 )
233*
234*              Apply the transformation as a rank-2 update:
235*                 A := A - v * w**T - w * v**T
236*
237               CALL DSPR2( UPLO, I, -ONE, AP( I1 ), 1, TAU, 1, AP )
238*
239               AP( I1+I-1 ) = E( I )
240            END IF
241            D( I+1 ) = AP( I1+I )
242            TAU( I ) = TAUI
243            I1 = I1 - I
244   10    CONTINUE
245         D( 1 ) = AP( 1 )
246      ELSE
247*
248*        Reduce the lower triangle of A. II is the index in AP of
249*        A(i,i) and I1I1 is the index of A(i+1,i+1).
250*
251         II = 1
252         DO 20 I = 1, N - 1
253            I1I1 = II + N - I + 1
254*
255*           Generate elementary reflector H(i) = I - tau * v * v**T
256*           to annihilate A(i+2:n,i)
257*
258            CALL DLARFG( N-I, AP( II+1 ), AP( II+2 ), 1, TAUI )
259            E( I ) = AP( II+1 )
260*
261            IF( TAUI.NE.ZERO ) THEN
262*
263*              Apply H(i) from both sides to A(i+1:n,i+1:n)
264*
265               AP( II+1 ) = ONE
266*
267*              Compute  y := tau * A * v  storing y in TAU(i:n-1)
268*
269               CALL DSPMV( UPLO, N-I, TAUI, AP( I1I1 ), AP( II+1 ), 1,
270     $                     ZERO, TAU( I ), 1 )
271*
272*              Compute  w := y - 1/2 * tau * (y**T *v) * v
273*
274               ALPHA = -HALF*TAUI*DDOT( N-I, TAU( I ), 1, AP( II+1 ),
275     $                 1 )
276               CALL DAXPY( N-I, ALPHA, AP( II+1 ), 1, TAU( I ), 1 )
277*
278*              Apply the transformation as a rank-2 update:
279*                 A := A - v * w**T - w * v**T
280*
281               CALL DSPR2( UPLO, N-I, -ONE, AP( II+1 ), 1, TAU( I ), 1,
282     $                     AP( I1I1 ) )
283*
284               AP( II+1 ) = E( I )
285            END IF
286            D( I ) = AP( II )
287            TAU( I ) = TAUI
288            II = I1I1
289   20    CONTINUE
290         D( N ) = AP( II )
291      END IF
292*
293      RETURN
294*
295*     End of DSPTRD
296*
297      END
298