1: *> \brief \b DGGSVP
2: *
3: * =========== DOCUMENTATION ===========
4: *
5: * Online html documentation available at
6: * http://www.netlib.org/lapack/explore-html/
7: *
8: *> \htmlonly
9: *> Download DGGSVP + dependencies
10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dggsvp.f">
11: *> [TGZ]</a>
12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dggsvp.f">
13: *> [ZIP]</a>
14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dggsvp.f">
15: *> [TXT]</a>
16: *> \endhtmlonly
17: *
18: * Definition:
19: * ===========
20: *
21: * SUBROUTINE DGGSVP( JOBU, JOBV, JOBQ, M, P, N, A, LDA, B, LDB,
22: * TOLA, TOLB, K, L, U, LDU, V, LDV, Q, LDQ,
23: * IWORK, TAU, WORK, INFO )
24: *
25: * .. Scalar Arguments ..
26: * CHARACTER JOBQ, JOBU, JOBV
27: * INTEGER INFO, K, L, LDA, LDB, LDQ, LDU, LDV, M, N, P
28: * DOUBLE PRECISION TOLA, TOLB
29: * ..
30: * .. Array Arguments ..
31: * INTEGER IWORK( * )
32: * DOUBLE PRECISION A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
33: * $ TAU( * ), U( LDU, * ), V( LDV, * ), WORK( * )
34: * ..
35: *
36: *
37: *> \par Purpose:
38: * =============
39: *>
40: *> \verbatim
41: *>
42: *> This routine is deprecated and has been replaced by routine DGGSVP3.
43: *>
44: *> DGGSVP computes orthogonal matrices U, V and Q such that
45: *>
46: *> N-K-L K L
47: *> U**T*A*Q = K ( 0 A12 A13 ) if M-K-L >= 0;
48: *> L ( 0 0 A23 )
49: *> M-K-L ( 0 0 0 )
50: *>
51: *> N-K-L K L
52: *> = K ( 0 A12 A13 ) if M-K-L < 0;
53: *> M-K ( 0 0 A23 )
54: *>
55: *> N-K-L K L
56: *> V**T*B*Q = L ( 0 0 B13 )
57: *> P-L ( 0 0 0 )
58: *>
59: *> where the K-by-K matrix A12 and L-by-L matrix B13 are nonsingular
60: *> upper triangular; A23 is L-by-L upper triangular if M-K-L >= 0,
61: *> otherwise A23 is (M-K)-by-L upper trapezoidal. K+L = the effective
62: *> numerical rank of the (M+P)-by-N matrix (A**T,B**T)**T.
63: *>
64: *> This decomposition is the preprocessing step for computing the
65: *> Generalized Singular Value Decomposition (GSVD), see subroutine
66: *> DGGSVD.
67: *> \endverbatim
68: *
69: * Arguments:
70: * ==========
71: *
72: *> \param[in] JOBU
73: *> \verbatim
74: *> JOBU is CHARACTER*1
75: *> = 'U': Orthogonal matrix U is computed;
76: *> = 'N': U is not computed.
77: *> \endverbatim
78: *>
79: *> \param[in] JOBV
80: *> \verbatim
81: *> JOBV is CHARACTER*1
82: *> = 'V': Orthogonal matrix V is computed;
83: *> = 'N': V is not computed.
84: *> \endverbatim
85: *>
86: *> \param[in] JOBQ
87: *> \verbatim
88: *> JOBQ is CHARACTER*1
89: *> = 'Q': Orthogonal matrix Q is computed;
90: *> = 'N': Q is not computed.
91: *> \endverbatim
92: *>
93: *> \param[in] M
94: *> \verbatim
95: *> M is INTEGER
96: *> The number of rows of the matrix A. M >= 0.
97: *> \endverbatim
98: *>
99: *> \param[in] P
100: *> \verbatim
101: *> P is INTEGER
102: *> The number of rows of the matrix B. P >= 0.
103: *> \endverbatim
104: *>
105: *> \param[in] N
106: *> \verbatim
107: *> N is INTEGER
108: *> The number of columns of the matrices A and B. N >= 0.
109: *> \endverbatim
110: *>
111: *> \param[in,out] A
112: *> \verbatim
113: *> A is DOUBLE PRECISION array, dimension (LDA,N)
114: *> On entry, the M-by-N matrix A.
115: *> On exit, A contains the triangular (or trapezoidal) matrix
116: *> described in the Purpose section.
117: *> \endverbatim
118: *>
119: *> \param[in] LDA
120: *> \verbatim
121: *> LDA is INTEGER
122: *> The leading dimension of the array A. LDA >= max(1,M).
123: *> \endverbatim
124: *>
125: *> \param[in,out] B
126: *> \verbatim
127: *> B is DOUBLE PRECISION array, dimension (LDB,N)
128: *> On entry, the P-by-N matrix B.
129: *> On exit, B contains the triangular matrix described in
130: *> the Purpose section.
131: *> \endverbatim
132: *>
133: *> \param[in] LDB
134: *> \verbatim
135: *> LDB is INTEGER
136: *> The leading dimension of the array B. LDB >= max(1,P).
137: *> \endverbatim
138: *>
139: *> \param[in] TOLA
140: *> \verbatim
141: *> TOLA is DOUBLE PRECISION
142: *> \endverbatim
143: *>
144: *> \param[in] TOLB
145: *> \verbatim
146: *> TOLB is DOUBLE PRECISION
147: *>
148: *> TOLA and TOLB are the thresholds to determine the effective
149: *> numerical rank of matrix B and a subblock of A. Generally,
150: *> they are set to
151: *> TOLA = MAX(M,N)*norm(A)*MACHEPS,
152: *> TOLB = MAX(P,N)*norm(B)*MACHEPS.
153: *> The size of TOLA and TOLB may affect the size of backward
154: *> errors of the decomposition.
155: *> \endverbatim
156: *>
157: *> \param[out] K
158: *> \verbatim
159: *> K is INTEGER
160: *> \endverbatim
161: *>
162: *> \param[out] L
163: *> \verbatim
164: *> L is INTEGER
165: *>
166: *> On exit, K and L specify the dimension of the subblocks
167: *> described in Purpose section.
168: *> K + L = effective numerical rank of (A**T,B**T)**T.
169: *> \endverbatim
170: *>
171: *> \param[out] U
172: *> \verbatim
173: *> U is DOUBLE PRECISION array, dimension (LDU,M)
174: *> If JOBU = 'U', U contains the orthogonal matrix U.
175: *> If JOBU = 'N', U is not referenced.
176: *> \endverbatim
177: *>
178: *> \param[in] LDU
179: *> \verbatim
180: *> LDU is INTEGER
181: *> The leading dimension of the array U. LDU >= max(1,M) if
182: *> JOBU = 'U'; LDU >= 1 otherwise.
183: *> \endverbatim
184: *>
185: *> \param[out] V
186: *> \verbatim
187: *> V is DOUBLE PRECISION array, dimension (LDV,P)
188: *> If JOBV = 'V', V contains the orthogonal matrix V.
189: *> If JOBV = 'N', V is not referenced.
190: *> \endverbatim
191: *>
192: *> \param[in] LDV
193: *> \verbatim
194: *> LDV is INTEGER
195: *> The leading dimension of the array V. LDV >= max(1,P) if
196: *> JOBV = 'V'; LDV >= 1 otherwise.
197: *> \endverbatim
198: *>
199: *> \param[out] Q
200: *> \verbatim
201: *> Q is DOUBLE PRECISION array, dimension (LDQ,N)
202: *> If JOBQ = 'Q', Q contains the orthogonal matrix Q.
203: *> If JOBQ = 'N', Q is not referenced.
204: *> \endverbatim
205: *>
206: *> \param[in] LDQ
207: *> \verbatim
208: *> LDQ is INTEGER
209: *> The leading dimension of the array Q. LDQ >= max(1,N) if
210: *> JOBQ = 'Q'; LDQ >= 1 otherwise.
211: *> \endverbatim
212: *>
213: *> \param[out] IWORK
214: *> \verbatim
215: *> IWORK is INTEGER array, dimension (N)
216: *> \endverbatim
217: *>
218: *> \param[out] TAU
219: *> \verbatim
220: *> TAU is DOUBLE PRECISION array, dimension (N)
221: *> \endverbatim
222: *>
223: *> \param[out] WORK
224: *> \verbatim
225: *> WORK is DOUBLE PRECISION array, dimension (max(3*N,M,P))
226: *> \endverbatim
227: *>
228: *> \param[out] INFO
229: *> \verbatim
230: *> INFO is INTEGER
231: *> = 0: successful exit
232: *> < 0: if INFO = -i, the i-th argument had an illegal value.
233: *> \endverbatim
234: *
235: * Authors:
236: * ========
237: *
238: *> \author Univ. of Tennessee
239: *> \author Univ. of California Berkeley
240: *> \author Univ. of Colorado Denver
241: *> \author NAG Ltd.
242: *
243: *> \ingroup doubleOTHERcomputational
244: *
245: *> \par Further Details:
246: * =====================
247: *>
248: *> The subroutine uses LAPACK subroutine DGEQPF for the QR factorization
249: *> with column pivoting to detect the effective numerical rank of the
250: *> a matrix. It may be replaced by a better rank determination strategy.
251: *>
252: * =====================================================================
253: SUBROUTINE DGGSVP( JOBU, JOBV, JOBQ, M, P, N, A, LDA, B, LDB,
254: $ TOLA, TOLB, K, L, U, LDU, V, LDV, Q, LDQ,
255: $ IWORK, TAU, WORK, INFO )
256: *
257: * -- LAPACK computational routine --
258: * -- LAPACK is a software package provided by Univ. of Tennessee, --
259: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
260: *
261: * .. Scalar Arguments ..
262: CHARACTER JOBQ, JOBU, JOBV
263: INTEGER INFO, K, L, LDA, LDB, LDQ, LDU, LDV, M, N, P
264: DOUBLE PRECISION TOLA, TOLB
265: * ..
266: * .. Array Arguments ..
267: INTEGER IWORK( * )
268: DOUBLE PRECISION A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
269: $ TAU( * ), U( LDU, * ), V( LDV, * ), WORK( * )
270: * ..
271: *
272: * =====================================================================
273: *
274: * .. Parameters ..
275: DOUBLE PRECISION ZERO, ONE
276: PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
277: * ..
278: * .. Local Scalars ..
279: LOGICAL FORWRD, WANTQ, WANTU, WANTV
280: INTEGER I, J
281: * ..
282: * .. External Functions ..
283: LOGICAL LSAME
284: EXTERNAL LSAME
285: * ..
286: * .. External Subroutines ..
287: EXTERNAL DGEQPF, DGEQR2, DGERQ2, DLACPY, DLAPMT, DLASET,
288: $ DORG2R, DORM2R, DORMR2, XERBLA
289: * ..
290: * .. Intrinsic Functions ..
291: INTRINSIC ABS, MAX, MIN
292: * ..
293: * .. Executable Statements ..
294: *
295: * Test the input parameters
296: *
297: WANTU = LSAME( JOBU, 'U' )
298: WANTV = LSAME( JOBV, 'V' )
299: WANTQ = LSAME( JOBQ, 'Q' )
300: FORWRD = .TRUE.
301: *
302: INFO = 0
303: IF( .NOT.( WANTU .OR. LSAME( JOBU, 'N' ) ) ) THEN
304: INFO = -1
305: ELSE IF( .NOT.( WANTV .OR. LSAME( JOBV, 'N' ) ) ) THEN
306: INFO = -2
307: ELSE IF( .NOT.( WANTQ .OR. LSAME( JOBQ, 'N' ) ) ) THEN
308: INFO = -3
309: ELSE IF( M.LT.0 ) THEN
310: INFO = -4
311: ELSE IF( P.LT.0 ) THEN
312: INFO = -5
313: ELSE IF( N.LT.0 ) THEN
314: INFO = -6
315: ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
316: INFO = -8
317: ELSE IF( LDB.LT.MAX( 1, P ) ) THEN
318: INFO = -10
319: ELSE IF( LDU.LT.1 .OR. ( WANTU .AND. LDU.LT.M ) ) THEN
320: INFO = -16
321: ELSE IF( LDV.LT.1 .OR. ( WANTV .AND. LDV.LT.P ) ) THEN
322: INFO = -18
323: ELSE IF( LDQ.LT.1 .OR. ( WANTQ .AND. LDQ.LT.N ) ) THEN
324: INFO = -20
325: END IF
326: IF( INFO.NE.0 ) THEN
327: CALL XERBLA( 'DGGSVP', -INFO )
328: RETURN
329: END IF
330: *
331: * QR with column pivoting of B: B*P = V*( S11 S12 )
332: * ( 0 0 )
333: *
334: DO 10 I = 1, N
335: IWORK( I ) = 0
336: 10 CONTINUE
337: CALL DGEQPF( P, N, B, LDB, IWORK, TAU, WORK, INFO )
338: *
339: * Update A := A*P
340: *
341: CALL DLAPMT( FORWRD, M, N, A, LDA, IWORK )
342: *
343: * Determine the effective rank of matrix B.
344: *
345: L = 0
346: DO 20 I = 1, MIN( P, N )
347: IF( ABS( B( I, I ) ).GT.TOLB )
348: $ L = L + 1
349: 20 CONTINUE
350: *
351: IF( WANTV ) THEN
352: *
353: * Copy the details of V, and form V.
354: *
355: CALL DLASET( 'Full', P, P, ZERO, ZERO, V, LDV )
356: IF( P.GT.1 )
357: $ CALL DLACPY( 'Lower', P-1, N, B( 2, 1 ), LDB, V( 2, 1 ),
358: $ LDV )
359: CALL DORG2R( P, P, MIN( P, N ), V, LDV, TAU, WORK, INFO )
360: END IF
361: *
362: * Clean up B
363: *
364: DO 40 J = 1, L - 1
365: DO 30 I = J + 1, L
366: B( I, J ) = ZERO
367: 30 CONTINUE
368: 40 CONTINUE
369: IF( P.GT.L )
370: $ CALL DLASET( 'Full', P-L, N, ZERO, ZERO, B( L+1, 1 ), LDB )
371: *
372: IF( WANTQ ) THEN
373: *
374: * Set Q = I and Update Q := Q*P
375: *
376: CALL DLASET( 'Full', N, N, ZERO, ONE, Q, LDQ )
377: CALL DLAPMT( FORWRD, N, N, Q, LDQ, IWORK )
378: END IF
379: *
380: IF( P.GE.L .AND. N.NE.L ) THEN
381: *
382: * RQ factorization of (S11 S12): ( S11 S12 ) = ( 0 S12 )*Z
383: *
384: CALL DGERQ2( L, N, B, LDB, TAU, WORK, INFO )
385: *
386: * Update A := A*Z**T
387: *
388: CALL DORMR2( 'Right', 'Transpose', M, N, L, B, LDB, TAU, A,
389: $ LDA, WORK, INFO )
390: *
391: IF( WANTQ ) THEN
392: *
393: * Update Q := Q*Z**T
394: *
395: CALL DORMR2( 'Right', 'Transpose', N, N, L, B, LDB, TAU, Q,
396: $ LDQ, WORK, INFO )
397: END IF
398: *
399: * Clean up B
400: *
401: CALL DLASET( 'Full', L, N-L, ZERO, ZERO, B, LDB )
402: DO 60 J = N - L + 1, N
403: DO 50 I = J - N + L + 1, L
404: B( I, J ) = ZERO
405: 50 CONTINUE
406: 60 CONTINUE
407: *
408: END IF
409: *
410: * Let N-L L
411: * A = ( A11 A12 ) M,
412: *
413: * then the following does the complete QR decomposition of A11:
414: *
415: * A11 = U*( 0 T12 )*P1**T
416: * ( 0 0 )
417: *
418: DO 70 I = 1, N - L
419: IWORK( I ) = 0
420: 70 CONTINUE
421: CALL DGEQPF( M, N-L, A, LDA, IWORK, TAU, WORK, INFO )
422: *
423: * Determine the effective rank of A11
424: *
425: K = 0
426: DO 80 I = 1, MIN( M, N-L )
427: IF( ABS( A( I, I ) ).GT.TOLA )
428: $ K = K + 1
429: 80 CONTINUE
430: *
431: * Update A12 := U**T*A12, where A12 = A( 1:M, N-L+1:N )
432: *
433: CALL DORM2R( 'Left', 'Transpose', M, L, MIN( M, N-L ), A, LDA,
434: $ TAU, A( 1, N-L+1 ), LDA, WORK, INFO )
435: *
436: IF( WANTU ) THEN
437: *
438: * Copy the details of U, and form U
439: *
440: CALL DLASET( 'Full', M, M, ZERO, ZERO, U, LDU )
441: IF( M.GT.1 )
442: $ CALL DLACPY( 'Lower', M-1, N-L, A( 2, 1 ), LDA, U( 2, 1 ),
443: $ LDU )
444: CALL DORG2R( M, M, MIN( M, N-L ), U, LDU, TAU, WORK, INFO )
445: END IF
446: *
447: IF( WANTQ ) THEN
448: *
449: * Update Q( 1:N, 1:N-L ) = Q( 1:N, 1:N-L )*P1
450: *
451: CALL DLAPMT( FORWRD, N, N-L, Q, LDQ, IWORK )
452: END IF
453: *
454: * Clean up A: set the strictly lower triangular part of
455: * A(1:K, 1:K) = 0, and A( K+1:M, 1:N-L ) = 0.
456: *
457: DO 100 J = 1, K - 1
458: DO 90 I = J + 1, K
459: A( I, J ) = ZERO
460: 90 CONTINUE
461: 100 CONTINUE
462: IF( M.GT.K )
463: $ CALL DLASET( 'Full', M-K, N-L, ZERO, ZERO, A( K+1, 1 ), LDA )
464: *
465: IF( N-L.GT.K ) THEN
466: *
467: * RQ factorization of ( T11 T12 ) = ( 0 T12 )*Z1
468: *
469: CALL DGERQ2( K, N-L, A, LDA, TAU, WORK, INFO )
470: *
471: IF( WANTQ ) THEN
472: *
473: * Update Q( 1:N,1:N-L ) = Q( 1:N,1:N-L )*Z1**T
474: *
475: CALL DORMR2( 'Right', 'Transpose', N, N-L, K, A, LDA, TAU,
476: $ Q, LDQ, WORK, INFO )
477: END IF
478: *
479: * Clean up A
480: *
481: CALL DLASET( 'Full', K, N-L-K, ZERO, ZERO, A, LDA )
482: DO 120 J = N - L - K + 1, N - L
483: DO 110 I = J - N + L + K + 1, K
484: A( I, J ) = ZERO
485: 110 CONTINUE
486: 120 CONTINUE
487: *
488: END IF
489: *
490: IF( M.GT.K ) THEN
491: *
492: * QR factorization of A( K+1:M,N-L+1:N )
493: *
494: CALL DGEQR2( M-K, L, A( K+1, N-L+1 ), LDA, TAU, WORK, INFO )
495: *
496: IF( WANTU ) THEN
497: *
498: * Update U(:,K+1:M) := U(:,K+1:M)*U1
499: *
500: CALL DORM2R( 'Right', 'No transpose', M, M-K, MIN( M-K, L ),
501: $ A( K+1, N-L+1 ), LDA, TAU, U( 1, K+1 ), LDU,
502: $ WORK, INFO )
503: END IF
504: *
505: * Clean up
506: *
507: DO 140 J = N - L + 1, N
508: DO 130 I = J - N + K + L + 1, M
509: A( I, J ) = ZERO
510: 130 CONTINUE
511: 140 CONTINUE
512: *
513: END IF
514: *
515: RETURN
516: *
517: * End of DGGSVP
518: *
519: END
CVSweb interface <joel.bertrand@systella.fr>