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