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