1*> \brief \b CGEBD2 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 CGEBD2 + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cgebd2.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cgebd2.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cgebd2.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18*  Definition:
19*  ===========
20*
21*       SUBROUTINE CGEBD2( 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               D( * ), E( * )
28*       COMPLEX            A( LDA, * ), TAUP( * ), TAUQ( * ), WORK( * )
29*       ..
30*
31*
32*> \par Purpose:
33*  =============
34*>
35*> \verbatim
36*>
37*> CGEBD2 reduces a complex general m by n matrix A to upper or lower
38*> real bidiagonal form B by a unitary transformation: Q**H * 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 COMPLEX 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 unitary matrix Q as a product of elementary
67*>            reflectors, and the elements above the first superdiagonal,
68*>            with the array TAUP, represent the unitary 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 unitary matrix Q as a product of
74*>            elementary reflectors, and the elements above the diagonal,
75*>            with the array TAUP, represent the unitary 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 COMPLEX array, dimension (min(M,N))
104*>          The scalar factors of the elementary reflectors which
105*>          represent the unitary matrix Q. See Further Details.
106*> \endverbatim
107*>
108*> \param[out] TAUP
109*> \verbatim
110*>          TAUP is COMPLEX array, dimension (min(M,N))
111*>          The scalar factors of the elementary reflectors which
112*>          represent the unitary matrix P. See Further Details.
113*> \endverbatim
114*>
115*> \param[out] WORK
116*> \verbatim
117*>          WORK is COMPLEX 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 complexGEcomputational
136* @precisions normal c -> s d z
137*
138*> \par Further Details:
139*  =====================
140*>
141*> \verbatim
142*>
143*>  The matrices Q and P are represented as products of elementary
144*>  reflectors:
145*>
146*>  If m >= n,
147*>
148*>     Q = H(1) H(2) . . . H(n)  and  P = G(1) G(2) . . . G(n-1)
149*>
150*>  Each H(i) and G(i) has the form:
151*>
152*>     H(i) = I - tauq * v * v**H  and G(i) = I - taup * u * u**H
153*>
154*>  where tauq and taup are complex scalars, and v and u are complex
155*>  vectors; v(1:i-1) = 0, v(i) = 1, and v(i+1:m) is stored on exit in
156*>  A(i+1:m,i); u(1:i) = 0, u(i+1) = 1, and u(i+2:n) is stored on exit in
157*>  A(i,i+2:n); tauq is stored in TAUQ(i) and taup in TAUP(i).
158*>
159*>  If m < n,
160*>
161*>     Q = H(1) H(2) . . . H(m-1)  and  P = G(1) G(2) . . . G(m)
162*>
163*>  Each H(i) and G(i) has the form:
164*>
165*>     H(i) = I - tauq * v * v**H  and G(i) = I - taup * u * u**H
166*>
167*>  where tauq and taup are complex scalars, v and u are complex vectors;
168*>  v(1:i) = 0, v(i+1) = 1, and v(i+2:m) is stored on exit in A(i+2:m,i);
169*>  u(1:i-1) = 0, u(i) = 1, and u(i+1:n) is stored on exit in A(i,i+1:n);
170*>  tauq is stored in TAUQ(i) and taup in TAUP(i).
171*>
172*>  The contents of A on exit are illustrated by the following examples:
173*>
174*>  m = 6 and n = 5 (m > n):          m = 5 and n = 6 (m < n):
175*>
176*>    (  d   e   u1  u1  u1 )           (  d   u1  u1  u1  u1  u1 )
177*>    (  v1  d   e   u2  u2 )           (  e   d   u2  u2  u2  u2 )
178*>    (  v1  v2  d   e   u3 )           (  v1  e   d   u3  u3  u3 )
179*>    (  v1  v2  v3  d   e  )           (  v1  v2  e   d   u4  u4 )
180*>    (  v1  v2  v3  v4  d  )           (  v1  v2  v3  e   d   u5 )
181*>    (  v1  v2  v3  v4  v5 )
182*>
183*>  where d and e denote diagonal and off-diagonal elements of B, vi
184*>  denotes an element of the vector defining H(i), and ui an element of
185*>  the vector defining G(i).
186*> \endverbatim
187*>
188*  =====================================================================
189      SUBROUTINE CGEBD2( M, N, A, LDA, D, E, TAUQ, TAUP, WORK, INFO )
190*
191*  -- LAPACK computational routine --
192*  -- LAPACK is a software package provided by Univ. of Tennessee,    --
193*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
194*
195*     .. Scalar Arguments ..
196      INTEGER            INFO, LDA, M, N
197*     ..
198*     .. Array Arguments ..
199      REAL               D( * ), E( * )
200      COMPLEX            A( LDA, * ), TAUP( * ), TAUQ( * ), WORK( * )
201*     ..
202*
203*  =====================================================================
204*
205*     .. Parameters ..
206      COMPLEX            ZERO, ONE
207      PARAMETER          ( ZERO = ( 0.0E+0, 0.0E+0 ),
208     $                   ONE = ( 1.0E+0, 0.0E+0 ) )
209*     ..
210*     .. Local Scalars ..
211      INTEGER            I
212      COMPLEX            ALPHA
213*     ..
214*     .. External Subroutines ..
215      EXTERNAL           CLACGV, CLARF, CLARFG, XERBLA
216*     ..
217*     .. Intrinsic Functions ..
218      INTRINSIC          CONJG, MAX, MIN
219*     ..
220*     .. Executable Statements ..
221*
222*     Test the input parameters
223*
224      INFO = 0
225      IF( M.LT.0 ) THEN
226         INFO = -1
227      ELSE IF( N.LT.0 ) THEN
228         INFO = -2
229      ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
230         INFO = -4
231      END IF
232      IF( INFO.LT.0 ) THEN
233         CALL XERBLA( 'CGEBD2', -INFO )
234         RETURN
235      END IF
236*
237      IF( M.GE.N ) THEN
238*
239*        Reduce to upper bidiagonal form
240*
241         DO 10 I = 1, N
242*
243*           Generate elementary reflector H(i) to annihilate A(i+1:m,i)
244*
245            ALPHA = A( I, I )
246            CALL CLARFG( M-I+1, ALPHA, A( MIN( I+1, M ), I ), 1,
247     $                   TAUQ( I ) )
248            D( I ) = REAL( ALPHA )
249            A( I, I ) = ONE
250*
251*           Apply H(i)**H to A(i:m,i+1:n) from the left
252*
253            IF( I.LT.N )
254     $         CALL CLARF( 'Left', M-I+1, N-I, A( I, I ), 1,
255     $                     CONJG( TAUQ( I ) ), A( I, I+1 ), LDA, WORK )
256            A( I, I ) = D( I )
257*
258            IF( I.LT.N ) THEN
259*
260*              Generate elementary reflector G(i) to annihilate
261*              A(i,i+2:n)
262*
263               CALL CLACGV( N-I, A( I, I+1 ), LDA )
264               ALPHA = A( I, I+1 )
265               CALL CLARFG( N-I, ALPHA, A( I, MIN( I+2, N ) ),
266     $                      LDA, TAUP( I ) )
267               E( I ) = REAL( ALPHA )
268               A( I, I+1 ) = ONE
269*
270*              Apply G(i) to A(i+1:m,i+1:n) from the right
271*
272               CALL CLARF( 'Right', M-I, N-I, A( I, I+1 ), LDA,
273     $                     TAUP( I ), A( I+1, I+1 ), LDA, WORK )
274               CALL CLACGV( N-I, A( I, I+1 ), LDA )
275               A( I, I+1 ) = E( I )
276            ELSE
277               TAUP( I ) = ZERO
278            END IF
279   10    CONTINUE
280      ELSE
281*
282*        Reduce to lower bidiagonal form
283*
284         DO 20 I = 1, M
285*
286*           Generate elementary reflector G(i) to annihilate A(i,i+1:n)
287*
288            CALL CLACGV( N-I+1, A( I, I ), LDA )
289            ALPHA = A( I, I )
290            CALL CLARFG( N-I+1, ALPHA, A( I, MIN( I+1, N ) ), LDA,
291     $                   TAUP( I ) )
292            D( I ) = REAL( ALPHA )
293            A( I, I ) = ONE
294*
295*           Apply G(i) to A(i+1:m,i:n) from the right
296*
297            IF( I.LT.M )
298     $         CALL CLARF( 'Right', M-I, N-I+1, A( I, I ), LDA,
299     $                     TAUP( I ), A( I+1, I ), LDA, WORK )
300            CALL CLACGV( N-I+1, A( I, I ), LDA )
301            A( I, I ) = D( I )
302*
303            IF( I.LT.M ) THEN
304*
305*              Generate elementary reflector H(i) to annihilate
306*              A(i+2:m,i)
307*
308               ALPHA = A( I+1, I )
309               CALL CLARFG( M-I, ALPHA, A( MIN( I+2, M ), I ), 1,
310     $                      TAUQ( I ) )
311               E( I ) = REAL( ALPHA )
312               A( I+1, I ) = ONE
313*
314*              Apply H(i)**H to A(i+1:m,i+1:n) from the left
315*
316               CALL CLARF( 'Left', M-I, N-I, A( I+1, I ), 1,
317     $                     CONJG( TAUQ( I ) ), A( I+1, I+1 ), LDA,
318     $                     WORK )
319               A( I+1, I ) = E( I )
320            ELSE
321               TAUQ( I ) = ZERO
322            END IF
323   20    CONTINUE
324      END IF
325      RETURN
326*
327*     End of CGEBD2
328*
329      END
330