1: *> \brief \b DORBDB4
2: *
3: * =========== DOCUMENTATION ===========
4: *
5: * Online html documentation available at
6: * http://www.netlib.org/lapack/explore-html/
7: *
8: *> \htmlonly
9: *> Download DORBDB4 + dependencies
10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dorbdb4.f">
11: *> [TGZ]</a>
12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dorbdb4.f">
13: *> [ZIP]</a>
14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dorbdb4.f">
15: *> [TXT]</a>
16: *> \endhtmlonly
17: *
18: * Definition:
19: * ===========
20: *
21: * SUBROUTINE DORBDB4( 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: * DOUBLE PRECISION PHANTOM(*), TAUP1(*), TAUP2(*), TAUQ1(*),
31: * $ WORK(*), X11(LDX11,*), X21(LDX21,*)
32: * ..
33: *
34: *
35: *> \par Purpose:
36: * =============
37: *>
38: *>\verbatim
39: *>
40: *> DORBDB4 simultaneously bidiagonalizes the blocks of a tall and skinny
41: *> matrix X with orthonomal 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 DORBDB1, DORBDB2, and DORBDB3 handle cases in
51: *> which M-Q is not the minimum dimension.
52: *>
53: *> The orthogonal 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 doubleOTHERcomputational
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 DORCSD for details.
197: *>
198: *> P1, P2, and Q1 are represented as products of elementary reflectors.
199: *> See DORCSD2BY1 for details on generating P1, P2, and Q1 using DORGQR
200: *> and DORGLQ.
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 DORBDB4( 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: DOUBLE PRECISION PHANTOM(*), TAUP1(*), TAUP2(*), TAUQ1(*),
224: $ WORK(*), X11(LDX11,*), X21(LDX21,*)
225: * ..
226: *
227: * ====================================================================
228: *
229: * .. Parameters ..
230: DOUBLE PRECISION NEGONE, ONE, ZERO
231: PARAMETER ( NEGONE = -1.0D0, ONE = 1.0D0, ZERO = 0.0D0 )
232: * ..
233: * .. Local Scalars ..
234: DOUBLE PRECISION C, S
235: INTEGER CHILDINFO, I, ILARF, IORBDB5, J, LLARF,
236: $ LORBDB5, LWORKMIN, LWORKOPT
237: LOGICAL LQUERY
238: * ..
239: * .. External Subroutines ..
240: EXTERNAL DLARF, DLARFGP, DORBDB5, DROT, DSCAL, XERBLA
241: * ..
242: * .. External Functions ..
243: DOUBLE PRECISION DNRM2
244: EXTERNAL DNRM2
245: * ..
246: * .. Intrinsic Function ..
247: INTRINSIC ATAN2, COS, MAX, SIN, SQRT
248: * ..
249: * .. Executable Statements ..
250: *
251: * Test input arguments
252: *
253: INFO = 0
254: LQUERY = LWORK .EQ. -1
255: *
256: IF( M .LT. 0 ) THEN
257: INFO = -1
258: ELSE IF( P .LT. M-Q .OR. M-P .LT. M-Q ) THEN
259: INFO = -2
260: ELSE IF( Q .LT. M-Q .OR. Q .GT. M ) THEN
261: INFO = -3
262: ELSE IF( LDX11 .LT. MAX( 1, P ) ) THEN
263: INFO = -5
264: ELSE IF( LDX21 .LT. MAX( 1, M-P ) ) THEN
265: INFO = -7
266: END IF
267: *
268: * Compute workspace
269: *
270: IF( INFO .EQ. 0 ) THEN
271: ILARF = 2
272: LLARF = MAX( Q-1, P-1, M-P-1 )
273: IORBDB5 = 2
274: LORBDB5 = Q
275: LWORKOPT = ILARF + LLARF - 1
276: LWORKOPT = MAX( LWORKOPT, IORBDB5 + LORBDB5 - 1 )
277: LWORKMIN = LWORKOPT
278: WORK(1) = LWORKOPT
279: IF( LWORK .LT. LWORKMIN .AND. .NOT.LQUERY ) THEN
280: INFO = -14
281: END IF
282: END IF
283: IF( INFO .NE. 0 ) THEN
284: CALL XERBLA( 'DORBDB4', -INFO )
285: RETURN
286: ELSE IF( LQUERY ) THEN
287: RETURN
288: END IF
289: *
290: * Reduce columns 1, ..., M-Q of X11 and X21
291: *
292: DO I = 1, M-Q
293: *
294: IF( I .EQ. 1 ) THEN
295: DO J = 1, M
296: PHANTOM(J) = ZERO
297: END DO
298: CALL DORBDB5( P, M-P, Q, PHANTOM(1), 1, PHANTOM(P+1), 1,
299: $ X11, LDX11, X21, LDX21, WORK(IORBDB5),
300: $ LORBDB5, CHILDINFO )
301: CALL DSCAL( P, NEGONE, PHANTOM(1), 1 )
302: CALL DLARFGP( P, PHANTOM(1), PHANTOM(2), 1, TAUP1(1) )
303: CALL DLARFGP( M-P, PHANTOM(P+1), PHANTOM(P+2), 1, TAUP2(1) )
304: THETA(I) = ATAN2( PHANTOM(1), PHANTOM(P+1) )
305: C = COS( THETA(I) )
306: S = SIN( THETA(I) )
307: PHANTOM(1) = ONE
308: PHANTOM(P+1) = ONE
309: CALL DLARF( 'L', P, Q, PHANTOM(1), 1, TAUP1(1), X11, LDX11,
310: $ WORK(ILARF) )
311: CALL DLARF( 'L', M-P, Q, PHANTOM(P+1), 1, TAUP2(1), X21,
312: $ LDX21, WORK(ILARF) )
313: ELSE
314: CALL DORBDB5( P-I+1, M-P-I+1, Q-I+1, X11(I,I-1), 1,
315: $ X21(I,I-1), 1, X11(I,I), LDX11, X21(I,I),
316: $ LDX21, WORK(IORBDB5), LORBDB5, CHILDINFO )
317: CALL DSCAL( P-I+1, NEGONE, X11(I,I-1), 1 )
318: CALL DLARFGP( P-I+1, X11(I,I-1), X11(I+1,I-1), 1, TAUP1(I) )
319: CALL DLARFGP( M-P-I+1, X21(I,I-1), X21(I+1,I-1), 1,
320: $ TAUP2(I) )
321: THETA(I) = ATAN2( X11(I,I-1), X21(I,I-1) )
322: C = COS( THETA(I) )
323: S = SIN( THETA(I) )
324: X11(I,I-1) = ONE
325: X21(I,I-1) = ONE
326: CALL DLARF( 'L', P-I+1, Q-I+1, X11(I,I-1), 1, TAUP1(I),
327: $ X11(I,I), LDX11, WORK(ILARF) )
328: CALL DLARF( 'L', M-P-I+1, Q-I+1, X21(I,I-1), 1, TAUP2(I),
329: $ X21(I,I), LDX21, WORK(ILARF) )
330: END IF
331: *
332: CALL DROT( Q-I+1, X11(I,I), LDX11, X21(I,I), LDX21, S, -C )
333: CALL DLARFGP( Q-I+1, X21(I,I), X21(I,I+1), LDX21, TAUQ1(I) )
334: C = X21(I,I)
335: X21(I,I) = ONE
336: CALL DLARF( 'R', P-I, Q-I+1, X21(I,I), LDX21, TAUQ1(I),
337: $ X11(I+1,I), LDX11, WORK(ILARF) )
338: CALL DLARF( 'R', M-P-I, Q-I+1, X21(I,I), LDX21, TAUQ1(I),
339: $ X21(I+1,I), LDX21, WORK(ILARF) )
340: IF( I .LT. M-Q ) THEN
341: S = SQRT( DNRM2( P-I, X11(I+1,I), 1 )**2
342: $ + DNRM2( M-P-I, X21(I+1,I), 1 )**2 )
343: PHI(I) = ATAN2( S, C )
344: END IF
345: *
346: END DO
347: *
348: * Reduce the bottom-right portion of X11 to [ I 0 ]
349: *
350: DO I = M - Q + 1, P
351: CALL DLARFGP( Q-I+1, X11(I,I), X11(I,I+1), LDX11, TAUQ1(I) )
352: X11(I,I) = ONE
353: CALL DLARF( 'R', P-I, Q-I+1, X11(I,I), LDX11, TAUQ1(I),
354: $ X11(I+1,I), LDX11, WORK(ILARF) )
355: CALL DLARF( 'R', Q-P, Q-I+1, X11(I,I), LDX11, TAUQ1(I),
356: $ X21(M-Q+1,I), LDX21, WORK(ILARF) )
357: END DO
358: *
359: * Reduce the bottom-right portion of X21 to [ 0 I ]
360: *
361: DO I = P + 1, Q
362: CALL DLARFGP( Q-I+1, X21(M-Q+I-P,I), X21(M-Q+I-P,I+1), LDX21,
363: $ TAUQ1(I) )
364: X21(M-Q+I-P,I) = ONE
365: CALL DLARF( 'R', Q-I, Q-I+1, X21(M-Q+I-P,I), LDX21, TAUQ1(I),
366: $ X21(M-Q+I-P+1,I), LDX21, WORK(ILARF) )
367: END DO
368: *
369: RETURN
370: *
371: * End of DORBDB4
372: *
373: END
374:
CVSweb interface <joel.bertrand@systella.fr>