LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cunbdb3.f
Go to the documentation of this file.
1*> \brief \b CUNBDB3
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download CUNBDB3 + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cunbdb3.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cunbdb3.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cunbdb3.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE CUNBDB3( M, P, Q, X11, LDX11, X21, LDX21, THETA, PHI,
22* TAUP1, TAUP2, TAUQ1, WORK, LWORK, INFO )
23*
24* .. Scalar Arguments ..
25* INTEGER INFO, LWORK, M, P, Q, LDX11, LDX21
26* ..
27* .. Array Arguments ..
28* REAL PHI(*), THETA(*)
29* COMPLEX TAUP1(*), TAUP2(*), TAUQ1(*), WORK(*),
30* $ X11(LDX11,*), X21(LDX21,*)
31* ..
32*
33*
34*> \par Purpose:
35* =============
36*>
37*>\verbatim
38*>
39*> CUNBDB3 simultaneously bidiagonalizes the blocks of a tall and skinny
40*> matrix X with orthonormal columns:
41*>
42*> [ B11 ]
43*> [ X11 ] [ P1 | ] [ 0 ]
44*> [-----] = [---------] [-----] Q1**T .
45*> [ X21 ] [ | P2 ] [ B21 ]
46*> [ 0 ]
47*>
48*> X11 is P-by-Q, and X21 is (M-P)-by-Q. M-P must be no larger than P,
49*> Q, or M-Q. Routines CUNBDB1, CUNBDB2, and CUNBDB4 handle cases in
50*> which M-P is not the minimum dimension.
51*>
52*> The unitary matrices P1, P2, and Q1 are P-by-P, (M-P)-by-(M-P),
53*> and (M-Q)-by-(M-Q), respectively. They are represented implicitly by
54*> Householder vectors.
55*>
56*> B11 and B12 are (M-P)-by-(M-P) bidiagonal matrices represented
57*> implicitly by angles THETA, PHI.
58*>
59*>\endverbatim
60*
61* Arguments:
62* ==========
63*
64*> \param[in] M
65*> \verbatim
66*> M is INTEGER
67*> The number of rows X11 plus the number of rows in X21.
68*> \endverbatim
69*>
70*> \param[in] P
71*> \verbatim
72*> P is INTEGER
73*> The number of rows in X11. 0 <= P <= M. M-P <= min(P,Q,M-Q).
74*> \endverbatim
75*>
76*> \param[in] Q
77*> \verbatim
78*> Q is INTEGER
79*> The number of columns in X11 and X21. 0 <= Q <= M.
80*> \endverbatim
81*>
82*> \param[in,out] X11
83*> \verbatim
84*> X11 is COMPLEX array, dimension (LDX11,Q)
85*> On entry, the top block of the matrix X to be reduced. On
86*> exit, the columns of tril(X11) specify reflectors for P1 and
87*> the rows of triu(X11,1) specify reflectors for Q1.
88*> \endverbatim
89*>
90*> \param[in] LDX11
91*> \verbatim
92*> LDX11 is INTEGER
93*> The leading dimension of X11. LDX11 >= P.
94*> \endverbatim
95*>
96*> \param[in,out] X21
97*> \verbatim
98*> X21 is COMPLEX array, dimension (LDX21,Q)
99*> On entry, the bottom block of the matrix X to be reduced. On
100*> exit, the columns of tril(X21) specify reflectors for P2.
101*> \endverbatim
102*>
103*> \param[in] LDX21
104*> \verbatim
105*> LDX21 is INTEGER
106*> The leading dimension of X21. LDX21 >= M-P.
107*> \endverbatim
108*>
109*> \param[out] THETA
110*> \verbatim
111*> THETA is REAL array, dimension (Q)
112*> The entries of the bidiagonal blocks B11, B21 are defined by
113*> THETA and PHI. See Further Details.
114*> \endverbatim
115*>
116*> \param[out] PHI
117*> \verbatim
118*> PHI is REAL array, dimension (Q-1)
119*> The entries of the bidiagonal blocks B11, B21 are defined by
120*> THETA and PHI. See Further Details.
121*> \endverbatim
122*>
123*> \param[out] TAUP1
124*> \verbatim
125*> TAUP1 is COMPLEX array, dimension (P)
126*> The scalar factors of the elementary reflectors that define
127*> P1.
128*> \endverbatim
129*>
130*> \param[out] TAUP2
131*> \verbatim
132*> TAUP2 is COMPLEX array, dimension (M-P)
133*> The scalar factors of the elementary reflectors that define
134*> P2.
135*> \endverbatim
136*>
137*> \param[out] TAUQ1
138*> \verbatim
139*> TAUQ1 is COMPLEX array, dimension (Q)
140*> The scalar factors of the elementary reflectors that define
141*> Q1.
142*> \endverbatim
143*>
144*> \param[out] WORK
145*> \verbatim
146*> WORK is COMPLEX array, dimension (LWORK)
147*> \endverbatim
148*>
149*> \param[in] LWORK
150*> \verbatim
151*> LWORK is INTEGER
152*> The dimension of the array WORK. LWORK >= M-Q.
153*>
154*> If LWORK = -1, then a workspace query is assumed; the routine
155*> only calculates the optimal size of the WORK array, returns
156*> this value as the first entry of the WORK array, and no error
157*> message related to LWORK is issued by XERBLA.
158*> \endverbatim
159*>
160*> \param[out] INFO
161*> \verbatim
162*> INFO is INTEGER
163*> = 0: successful exit.
164*> < 0: if INFO = -i, the i-th argument had an illegal value.
165*> \endverbatim
166*>
167*
168* Authors:
169* ========
170*
171*> \author Univ. of Tennessee
172*> \author Univ. of California Berkeley
173*> \author Univ. of Colorado Denver
174*> \author NAG Ltd.
175*
176*> \ingroup unbdb3
177*
178*> \par Further Details:
179* =====================
180
181*> \verbatim
182*>
183*> The upper-bidiagonal blocks B11, B21 are represented implicitly by
184*> angles THETA(1), ..., THETA(Q) and PHI(1), ..., PHI(Q-1). Every entry
185*> in each bidiagonal band is a product of a sine or cosine of a THETA
186*> with a sine or cosine of a PHI. See [1] or CUNCSD for details.
187*>
188*> P1, P2, and Q1 are represented as products of elementary reflectors.
189*> See CUNCSD2BY1 for details on generating P1, P2, and Q1 using CUNGQR
190*> and CUNGLQ.
191*> \endverbatim
192*
193*> \par References:
194* ================
195*>
196*> [1] Brian D. Sutton. Computing the complete CS decomposition. Numer.
197*> Algorithms, 50(1):33-65, 2009.
198*>
199* =====================================================================
200 SUBROUTINE cunbdb3( M, P, Q, X11, LDX11, X21, LDX21, THETA, PHI,
201 $ TAUP1, TAUP2, TAUQ1, WORK, LWORK, INFO )
202*
203* -- LAPACK computational routine --
204* -- LAPACK is a software package provided by Univ. of Tennessee, --
205* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
206*
207* .. Scalar Arguments ..
208 INTEGER INFO, LWORK, M, P, Q, LDX11, LDX21
209* ..
210* .. Array Arguments ..
211 REAL PHI(*), THETA(*)
212 COMPLEX TAUP1(*), TAUP2(*), TAUQ1(*), WORK(*),
213 $ x11(ldx11,*), x21(ldx21,*)
214* ..
215*
216* ====================================================================
217*
218* .. Parameters ..
219 COMPLEX ONE
220 parameter( one = (1.0e0,0.0e0) )
221* ..
222* .. Local Scalars ..
223 REAL C, S
224 INTEGER CHILDINFO, I, ILARF, IORBDB5, LLARF, LORBDB5,
225 $ lworkmin, lworkopt
226 LOGICAL LQUERY
227* ..
228* .. External Subroutines ..
229 EXTERNAL clarf, clarfgp, cunbdb5, csrot, clacgv, xerbla
230* ..
231* .. External Functions ..
232 REAL SCNRM2, SROUNDUP_LWORK
233 EXTERNAL scnrm2, sroundup_lwork
234* ..
235* .. Intrinsic Function ..
236 INTRINSIC atan2, cos, max, sin, sqrt
237* ..
238* .. Executable Statements ..
239*
240* Test input arguments
241*
242 info = 0
243 lquery = lwork .EQ. -1
244*
245 IF( m .LT. 0 ) THEN
246 info = -1
247 ELSE IF( 2*p .LT. m .OR. p .GT. m ) THEN
248 info = -2
249 ELSE IF( q .LT. m-p .OR. m-q .LT. m-p ) THEN
250 info = -3
251 ELSE IF( ldx11 .LT. max( 1, p ) ) THEN
252 info = -5
253 ELSE IF( ldx21 .LT. max( 1, m-p ) ) THEN
254 info = -7
255 END IF
256*
257* Compute workspace
258*
259 IF( info .EQ. 0 ) THEN
260 ilarf = 2
261 llarf = max( p, m-p-1, q-1 )
262 iorbdb5 = 2
263 lorbdb5 = q-1
264 lworkopt = max( ilarf+llarf-1, iorbdb5+lorbdb5-1 )
265 lworkmin = lworkopt
266 work(1) = sroundup_lwork(lworkopt)
267 IF( lwork .LT. lworkmin .AND. .NOT.lquery ) THEN
268 info = -14
269 END IF
270 END IF
271 IF( info .NE. 0 ) THEN
272 CALL xerbla( 'CUNBDB3', -info )
273 RETURN
274 ELSE IF( lquery ) THEN
275 RETURN
276 END IF
277*
278* Reduce rows 1, ..., M-P of X11 and X21
279*
280 DO i = 1, m-p
281*
282 IF( i .GT. 1 ) THEN
283 CALL csrot( q-i+1, x11(i-1,i), ldx11, x21(i,i), ldx11, c,
284 $ s )
285 END IF
286*
287 CALL clacgv( q-i+1, x21(i,i), ldx21 )
288 CALL clarfgp( q-i+1, x21(i,i), x21(i,i+1), ldx21, tauq1(i) )
289 s = real( x21(i,i) )
290 x21(i,i) = one
291 CALL clarf( 'R', p-i+1, q-i+1, x21(i,i), ldx21, tauq1(i),
292 $ x11(i,i), ldx11, work(ilarf) )
293 CALL clarf( 'R', m-p-i, q-i+1, x21(i,i), ldx21, tauq1(i),
294 $ x21(i+1,i), ldx21, work(ilarf) )
295 CALL clacgv( q-i+1, x21(i,i), ldx21 )
296 c = sqrt( scnrm2( p-i+1, x11(i,i), 1 )**2
297 $ + scnrm2( m-p-i, x21(i+1,i), 1 )**2 )
298 theta(i) = atan2( s, c )
299*
300 CALL cunbdb5( p-i+1, m-p-i, q-i, x11(i,i), 1, x21(i+1,i), 1,
301 $ x11(i,i+1), ldx11, x21(i+1,i+1), ldx21,
302 $ work(iorbdb5), lorbdb5, childinfo )
303 CALL clarfgp( p-i+1, x11(i,i), x11(i+1,i), 1, taup1(i) )
304 IF( i .LT. m-p ) THEN
305 CALL clarfgp( m-p-i, x21(i+1,i), x21(i+2,i), 1, taup2(i) )
306 phi(i) = atan2( real( x21(i+1,i) ), real( x11(i,i) ) )
307 c = cos( phi(i) )
308 s = sin( phi(i) )
309 x21(i+1,i) = one
310 CALL clarf( 'L', m-p-i, q-i, x21(i+1,i), 1, conjg(taup2(i)),
311 $ x21(i+1,i+1), ldx21, work(ilarf) )
312 END IF
313 x11(i,i) = one
314 CALL clarf( 'L', p-i+1, q-i, x11(i,i), 1, conjg(taup1(i)),
315 $ x11(i,i+1), ldx11, work(ilarf) )
316*
317 END DO
318*
319* Reduce the bottom-right portion of X11 to the identity matrix
320*
321 DO i = m-p + 1, q
322 CALL clarfgp( p-i+1, x11(i,i), x11(i+1,i), 1, taup1(i) )
323 x11(i,i) = one
324 CALL clarf( 'L', p-i+1, q-i, x11(i,i), 1, conjg(taup1(i)),
325 $ x11(i,i+1), ldx11, work(ilarf) )
326 END DO
327*
328 RETURN
329*
330* End of CUNBDB3
331*
332 END
333
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine clacgv(n, x, incx)
CLACGV conjugates a complex vector.
Definition clacgv.f:74
subroutine clarf(side, m, n, v, incv, tau, c, ldc, work)
CLARF applies an elementary reflector to a general rectangular matrix.
Definition clarf.f:128
subroutine clarfgp(n, alpha, x, incx, tau)
CLARFGP generates an elementary reflector (Householder matrix) with non-negative beta.
Definition clarfgp.f:104
subroutine csrot(n, cx, incx, cy, incy, c, s)
CSROT
Definition csrot.f:98
subroutine cunbdb3(m, p, q, x11, ldx11, x21, ldx21, theta, phi, taup1, taup2, tauq1, work, lwork, info)
CUNBDB3
Definition cunbdb3.f:202
subroutine cunbdb5(m1, m2, n, x1, incx1, x2, incx2, q1, ldq1, q2, ldq2, work, lwork, info)
CUNBDB5
Definition cunbdb5.f:156