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