1*> \brief \b SGEBD2 reduces a general matrix to bidiagonal form using an unblocked algorithm.
2*
3*  =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6*            http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download SGEBD2 + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgebd2.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgebd2.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgebd2.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18*  Definition:
19*  ===========
20*
21*       SUBROUTINE SGEBD2( M, N, A, LDA, D, E, TAUQ, TAUP, WORK, INFO )
22*
23*       .. Scalar Arguments ..
24*       INTEGER            INFO, LDA, M, N
25*       ..
26*       .. Array Arguments ..
27*       REAL               A( LDA, * ), D( * ), E( * ), TAUP( * ),
28*      $                   TAUQ( * ), WORK( * )
29*       ..
30*
31*
32*> \par Purpose:
33*  =============
34*>
35*> \verbatim
36*>
37*> SGEBD2 reduces a real general m by n matrix A to upper or lower
38*> bidiagonal form B by an orthogonal transformation: Q**T * A * P = B.
39*>
40*> If m >= n, B is upper bidiagonal; if m < n, B is lower bidiagonal.
41*> \endverbatim
42*
43*  Arguments:
44*  ==========
45*
46*> \param[in] M
47*> \verbatim
48*>          M is INTEGER
49*>          The number of rows in the matrix A.  M >= 0.
50*> \endverbatim
51*>
52*> \param[in] N
53*> \verbatim
54*>          N is INTEGER
55*>          The number of columns in the matrix A.  N >= 0.
56*> \endverbatim
57*>
58*> \param[in,out] A
59*> \verbatim
60*>          A is REAL array, dimension (LDA,N)
61*>          On entry, the m by n general matrix to be reduced.
62*>          On exit,
63*>          if m >= n, the diagonal and the first superdiagonal are
64*>            overwritten with the upper bidiagonal matrix B; the
65*>            elements below the diagonal, with the array TAUQ, represent
66*>            the orthogonal matrix Q as a product of elementary
67*>            reflectors, and the elements above the first superdiagonal,
68*>            with the array TAUP, represent the orthogonal matrix P as
69*>            a product of elementary reflectors;
70*>          if m < n, the diagonal and the first subdiagonal are
71*>            overwritten with the lower bidiagonal matrix B; the
72*>            elements below the first subdiagonal, with the array TAUQ,
73*>            represent the orthogonal matrix Q as a product of
74*>            elementary reflectors, and the elements above the diagonal,
75*>            with the array TAUP, represent the orthogonal matrix P as
76*>            a product of elementary reflectors.
77*>          See Further Details.
78*> \endverbatim
79*>
80*> \param[in] LDA
81*> \verbatim
82*>          LDA is INTEGER
83*>          The leading dimension of the array A.  LDA >= max(1,M).
84*> \endverbatim
85*>
86*> \param[out] D
87*> \verbatim
88*>          D is REAL array, dimension (min(M,N))
89*>          The diagonal elements of the bidiagonal matrix B:
90*>          D(i) = A(i,i).
91*> \endverbatim
92*>
93*> \param[out] E
94*> \verbatim
95*>          E is REAL array, dimension (min(M,N)-1)
96*>          The off-diagonal elements of the bidiagonal matrix B:
97*>          if m >= n, E(i) = A(i,i+1) for i = 1,2,...,n-1;
98*>          if m < n, E(i) = A(i+1,i) for i = 1,2,...,m-1.
99*> \endverbatim
100*>
101*> \param[out] TAUQ
102*> \verbatim
103*>          TAUQ is REAL array, dimension (min(M,N))
104*>          The scalar factors of the elementary reflectors which
105*>          represent the orthogonal matrix Q. See Further Details.
106*> \endverbatim
107*>
108*> \param[out] TAUP
109*> \verbatim
110*>          TAUP is REAL array, dimension (min(M,N))
111*>          The scalar factors of the elementary reflectors which
112*>          represent the orthogonal matrix P. See Further Details.
113*> \endverbatim
114*>
115*> \param[out] WORK
116*> \verbatim
117*>          WORK is REAL array, dimension (max(M,N))
118*> \endverbatim
119*>
120*> \param[out] INFO
121*> \verbatim
122*>          INFO is INTEGER
123*>          = 0: successful exit.
124*>          < 0: if INFO = -i, the i-th argument had an illegal value.
125*> \endverbatim
126*
127*  Authors:
128*  ========
129*
130*> \author Univ. of Tennessee
131*> \author Univ. of California Berkeley
132*> \author Univ. of Colorado Denver
133*> \author NAG Ltd.
134*
135*> \ingroup realGEcomputational
136*
137*> \par Further Details:
138*  =====================
139*>
140*> \verbatim
141*>
142*>  The matrices Q and P are represented as products of elementary
143*>  reflectors:
144*>
145*>  If m >= n,
146*>
147*>     Q = H(1) H(2) . . . H(n)  and  P = G(1) G(2) . . . G(n-1)
148*>
149*>  Each H(i) and G(i) has the form:
150*>
151*>     H(i) = I - tauq * v * v**T  and G(i) = I - taup * u * u**T
152*>
153*>  where tauq and taup are real scalars, and v and u are real vectors;
154*>  v(1:i-1) = 0, v(i) = 1, and v(i+1:m) is stored on exit in A(i+1:m,i);
155*>  u(1:i) = 0, u(i+1) = 1, and u(i+2:n) is stored on exit in A(i,i+2:n);
156*>  tauq is stored in TAUQ(i) and taup in TAUP(i).
157*>
158*>  If m < n,
159*>
160*>     Q = H(1) H(2) . . . H(m-1)  and  P = G(1) G(2) . . . G(m)
161*>
162*>  Each H(i) and G(i) has the form:
163*>
164*>     H(i) = I - tauq * v * v**T  and G(i) = I - taup * u * u**T
165*>
166*>  where tauq and taup are real scalars, and v and u are real vectors;
167*>  v(1:i) = 0, v(i+1) = 1, and v(i+2:m) is stored on exit in A(i+2:m,i);
168*>  u(1:i-1) = 0, u(i) = 1, and u(i+1:n) is stored on exit in A(i,i+1:n);
169*>  tauq is stored in TAUQ(i) and taup in TAUP(i).
170*>
171*>  The contents of A on exit are illustrated by the following examples:
172*>
173*>  m = 6 and n = 5 (m > n):          m = 5 and n = 6 (m < n):
174*>
175*>    (  d   e   u1  u1  u1 )           (  d   u1  u1  u1  u1  u1 )
176*>    (  v1  d   e   u2  u2 )           (  e   d   u2  u2  u2  u2 )
177*>    (  v1  v2  d   e   u3 )           (  v1  e   d   u3  u3  u3 )
178*>    (  v1  v2  v3  d   e  )           (  v1  v2  e   d   u4  u4 )
179*>    (  v1  v2  v3  v4  d  )           (  v1  v2  v3  e   d   u5 )
180*>    (  v1  v2  v3  v4  v5 )
181*>
182*>  where d and e denote diagonal and off-diagonal elements of B, vi
183*>  denotes an element of the vector defining H(i), and ui an element of
184*>  the vector defining G(i).
185*> \endverbatim
186*>
187*  =====================================================================
188      SUBROUTINE SGEBD2( M, N, A, LDA, D, E, TAUQ, TAUP, WORK, INFO )
189*
190*  -- LAPACK computational routine --
191*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
192*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
193*
194*     .. Scalar Arguments ..
195      INTEGER            INFO, LDA, M, N
196*     ..
197*     .. Array Arguments ..
198      REAL               A( LDA, * ), D( * ), E( * ), TAUP( * ),
199     $                   TAUQ( * ), WORK( * )
200*     ..
201*
202*  =====================================================================
203*
204*     .. Parameters ..
205      REAL               ZERO, ONE
206      PARAMETER          ( ZERO = 0.0E+0, ONE = 1.0E+0 )
207*     ..
208*     .. Local Scalars ..
209      INTEGER            I
210*     ..
211*     .. External Subroutines ..
212      EXTERNAL           SLARF, SLARFG, XERBLA
213*     ..
214*     .. Intrinsic Functions ..
215      INTRINSIC          MAX, MIN
216*     ..
217*     .. Executable Statements ..
218*
219*     Test the input parameters
220*
221      INFO = 0
222      IF( M.LT.0 ) THEN
223         INFO = -1
224      ELSE IF( N.LT.0 ) THEN
225         INFO = -2
226      ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
227         INFO = -4
228      END IF
229      IF( INFO.LT.0 ) THEN
230         CALL XERBLA( 'SGEBD2', -INFO )
231         RETURN
232      END IF
233*
234      IF( M.GE.N ) THEN
235*
236*        Reduce to upper bidiagonal form
237*
238         DO 10 I = 1, N
239*
240*           Generate elementary reflector H(i) to annihilate A(i+1:m,i)
241*
242            CALL SLARFG( M-I+1, A( I, I ), A( MIN( I+1, M ), I ), 1,
243     $                   TAUQ( I ) )
244            D( I ) = A( I, I )
245            A( I, I ) = ONE
246*
247*           Apply H(i) to A(i:m,i+1:n) from the left
248*
249            IF( I.LT.N )
250     $         CALL SLARF( 'Left', M-I+1, N-I, A( I, I ), 1, TAUQ( I ),
251     $                     A( I, I+1 ), LDA, WORK )
252            A( I, I ) = D( I )
253*
254            IF( I.LT.N ) THEN
255*
256*              Generate elementary reflector G(i) to annihilate
257*              A(i,i+2:n)
258*
259               CALL SLARFG( N-I, A( I, I+1 ), A( I, MIN( I+2, N ) ),
260     $                      LDA, TAUP( I ) )
261               E( I ) = A( I, I+1 )
262               A( I, I+1 ) = ONE
263*
264*              Apply G(i) to A(i+1:m,i+1:n) from the right
265*
266               CALL SLARF( 'Right', M-I, N-I, A( I, I+1 ), LDA,
267     $                     TAUP( I ), A( I+1, I+1 ), LDA, WORK )
268               A( I, I+1 ) = E( I )
269            ELSE
270               TAUP( I ) = ZERO
271            END IF
272   10    CONTINUE
273      ELSE
274*
275*        Reduce to lower bidiagonal form
276*
277         DO 20 I = 1, M
278*
279*           Generate elementary reflector G(i) to annihilate A(i,i+1:n)
280*
281            CALL SLARFG( N-I+1, A( I, I ), A( I, MIN( I+1, N ) ), LDA,
282     $                   TAUP( I ) )
283            D( I ) = A( I, I )
284            A( I, I ) = ONE
285*
286*           Apply G(i) to A(i+1:m,i:n) from the right
287*
288            IF( I.LT.M )
289     $         CALL SLARF( 'Right', M-I, N-I+1, A( I, I ), LDA,
290     $                     TAUP( I ), A( I+1, I ), LDA, WORK )
291            A( I, I ) = D( I )
292*
293            IF( I.LT.M ) THEN
294*
295*              Generate elementary reflector H(i) to annihilate
296*              A(i+2:m,i)
297*
298               CALL SLARFG( M-I, A( I+1, I ), A( MIN( I+2, M ), I ), 1,
299     $                      TAUQ( I ) )
300               E( I ) = A( I+1, I )
301               A( I+1, I ) = ONE
302*
303*              Apply H(i) to A(i+1:m,i+1:n) from the left
304*
305               CALL SLARF( 'Left', M-I, N-I, A( I+1, I ), 1, TAUQ( I ),
306     $                     A( I+1, I+1 ), LDA, WORK )
307               A( I+1, I ) = E( I )
308            ELSE
309               TAUQ( I ) = ZERO
310            END IF
311   20    CONTINUE
312      END IF
313      RETURN
314*
315*     End of SGEBD2
316*
317      END
318