1: *> \brief <b> ZGGSVD3 computes the singular value decomposition (SVD) for OTHER matrices</b>
2: *
3: * =========== DOCUMENTATION ===========
4: *
5: * Online html documentation available at
6: * http://www.netlib.org/lapack/explore-html/
7: *
8: *> \htmlonly
9: *> Download ZGGSVD3 + dependencies
10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zggsvd3.f">
11: *> [TGZ]</a>
12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zggsvd3.f">
13: *> [ZIP]</a>
14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zggsvd3.f">
15: *> [TXT]</a>
16: *> \endhtmlonly
17: *
18: * Definition:
19: * ===========
20: *
21: * SUBROUTINE ZGGSVD3( JOBU, JOBV, JOBQ, M, N, P, K, L, A, LDA, B,
22: * LDB, ALPHA, BETA, U, LDU, V, LDV, Q, LDQ, WORK,
23: * LWORK, RWORK, IWORK, INFO )
24: *
25: * .. Scalar Arguments ..
26: * CHARACTER JOBQ, JOBU, JOBV
27: * INTEGER INFO, K, L, LDA, LDB, LDQ, LDU, LDV, M, N, P, LWORK
28: * ..
29: * .. Array Arguments ..
30: * INTEGER IWORK( * )
31: * DOUBLE PRECISION ALPHA( * ), BETA( * ), RWORK( * )
32: * COMPLEX*16 A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
33: * $ U( LDU, * ), V( LDV, * ), WORK( * )
34: * ..
35: *
36: *
37: *> \par Purpose:
38: * =============
39: *>
40: *> \verbatim
41: *>
42: *> ZGGSVD3 computes the generalized singular value decomposition (GSVD)
43: *> of an M-by-N complex matrix A and P-by-N complex matrix B:
44: *>
45: *> U**H*A*Q = D1*( 0 R ), V**H*B*Q = D2*( 0 R )
46: *>
47: *> where U, V and Q are unitary matrices.
48: *> Let K+L = the effective numerical rank of the
49: *> matrix (A**H,B**H)**H, then R is a (K+L)-by-(K+L) nonsingular upper
50: *> triangular matrix, D1 and D2 are M-by-(K+L) and P-by-(K+L) "diagonal"
51: *> matrices and of the following structures, respectively:
52: *>
53: *> If M-K-L >= 0,
54: *>
55: *> K L
56: *> D1 = K ( I 0 )
57: *> L ( 0 C )
58: *> M-K-L ( 0 0 )
59: *>
60: *> K L
61: *> D2 = L ( 0 S )
62: *> P-L ( 0 0 )
63: *>
64: *> N-K-L K L
65: *> ( 0 R ) = K ( 0 R11 R12 )
66: *> L ( 0 0 R22 )
67: *> where
68: *>
69: *> C = diag( ALPHA(K+1), ... , ALPHA(K+L) ),
70: *> S = diag( BETA(K+1), ... , BETA(K+L) ),
71: *> C**2 + S**2 = I.
72: *>
73: *> R is stored in A(1:K+L,N-K-L+1:N) on exit.
74: *>
75: *> If M-K-L < 0,
76: *>
77: *> K M-K K+L-M
78: *> D1 = K ( I 0 0 )
79: *> M-K ( 0 C 0 )
80: *>
81: *> K M-K K+L-M
82: *> D2 = M-K ( 0 S 0 )
83: *> K+L-M ( 0 0 I )
84: *> P-L ( 0 0 0 )
85: *>
86: *> N-K-L K M-K K+L-M
87: *> ( 0 R ) = K ( 0 R11 R12 R13 )
88: *> M-K ( 0 0 R22 R23 )
89: *> K+L-M ( 0 0 0 R33 )
90: *>
91: *> where
92: *>
93: *> C = diag( ALPHA(K+1), ... , ALPHA(M) ),
94: *> S = diag( BETA(K+1), ... , BETA(M) ),
95: *> C**2 + S**2 = I.
96: *>
97: *> (R11 R12 R13 ) is stored in A(1:M, N-K-L+1:N), and R33 is stored
98: *> ( 0 R22 R23 )
99: *> in B(M-K+1:L,N+M-K-L+1:N) on exit.
100: *>
101: *> The routine computes C, S, R, and optionally the unitary
102: *> transformation matrices U, V and Q.
103: *>
104: *> In particular, if B is an N-by-N nonsingular matrix, then the GSVD of
105: *> A and B implicitly gives the SVD of A*inv(B):
106: *> A*inv(B) = U*(D1*inv(D2))*V**H.
107: *> If ( A**H,B**H)**H has orthonormal columns, then the GSVD of A and B is also
108: *> equal to the CS decomposition of A and B. Furthermore, the GSVD can
109: *> be used to derive the solution of the eigenvalue problem:
110: *> A**H*A x = lambda* B**H*B x.
111: *> In some literature, the GSVD of A and B is presented in the form
112: *> U**H*A*X = ( 0 D1 ), V**H*B*X = ( 0 D2 )
113: *> where U and V are orthogonal and X is nonsingular, and D1 and D2 are
114: *> ``diagonal''. The former GSVD form can be converted to the latter
115: *> form by taking the nonsingular matrix X as
116: *>
117: *> X = Q*( I 0 )
118: *> ( 0 inv(R) )
119: *> \endverbatim
120: *
121: * Arguments:
122: * ==========
123: *
124: *> \param[in] JOBU
125: *> \verbatim
126: *> JOBU is CHARACTER*1
127: *> = 'U': Unitary matrix U is computed;
128: *> = 'N': U is not computed.
129: *> \endverbatim
130: *>
131: *> \param[in] JOBV
132: *> \verbatim
133: *> JOBV is CHARACTER*1
134: *> = 'V': Unitary matrix V is computed;
135: *> = 'N': V is not computed.
136: *> \endverbatim
137: *>
138: *> \param[in] JOBQ
139: *> \verbatim
140: *> JOBQ is CHARACTER*1
141: *> = 'Q': Unitary matrix Q is computed;
142: *> = 'N': Q is not computed.
143: *> \endverbatim
144: *>
145: *> \param[in] M
146: *> \verbatim
147: *> M is INTEGER
148: *> The number of rows of the matrix A. M >= 0.
149: *> \endverbatim
150: *>
151: *> \param[in] N
152: *> \verbatim
153: *> N is INTEGER
154: *> The number of columns of the matrices A and B. N >= 0.
155: *> \endverbatim
156: *>
157: *> \param[in] P
158: *> \verbatim
159: *> P is INTEGER
160: *> The number of rows of the matrix B. P >= 0.
161: *> \endverbatim
162: *>
163: *> \param[out] K
164: *> \verbatim
165: *> K is INTEGER
166: *> \endverbatim
167: *>
168: *> \param[out] L
169: *> \verbatim
170: *> L is INTEGER
171: *>
172: *> On exit, K and L specify the dimension of the subblocks
173: *> described in Purpose.
174: *> K + L = effective numerical rank of (A**H,B**H)**H.
175: *> \endverbatim
176: *>
177: *> \param[in,out] A
178: *> \verbatim
179: *> A is COMPLEX*16 array, dimension (LDA,N)
180: *> On entry, the M-by-N matrix A.
181: *> On exit, A contains the triangular matrix R, or part of R.
182: *> See Purpose for details.
183: *> \endverbatim
184: *>
185: *> \param[in] LDA
186: *> \verbatim
187: *> LDA is INTEGER
188: *> The leading dimension of the array A. LDA >= max(1,M).
189: *> \endverbatim
190: *>
191: *> \param[in,out] B
192: *> \verbatim
193: *> B is COMPLEX*16 array, dimension (LDB,N)
194: *> On entry, the P-by-N matrix B.
195: *> On exit, B contains part of the triangular matrix R if
196: *> M-K-L < 0. See Purpose for details.
197: *> \endverbatim
198: *>
199: *> \param[in] LDB
200: *> \verbatim
201: *> LDB is INTEGER
202: *> The leading dimension of the array B. LDB >= max(1,P).
203: *> \endverbatim
204: *>
205: *> \param[out] ALPHA
206: *> \verbatim
207: *> ALPHA is DOUBLE PRECISION array, dimension (N)
208: *> \endverbatim
209: *>
210: *> \param[out] BETA
211: *> \verbatim
212: *> BETA is DOUBLE PRECISION array, dimension (N)
213: *>
214: *> On exit, ALPHA and BETA contain the generalized singular
215: *> value pairs of A and B;
216: *> ALPHA(1:K) = 1,
217: *> BETA(1:K) = 0,
218: *> and if M-K-L >= 0,
219: *> ALPHA(K+1:K+L) = C,
220: *> BETA(K+1:K+L) = S,
221: *> or if M-K-L < 0,
222: *> ALPHA(K+1:M)=C, ALPHA(M+1:K+L)=0
223: *> BETA(K+1:M) =S, BETA(M+1:K+L) =1
224: *> and
225: *> ALPHA(K+L+1:N) = 0
226: *> BETA(K+L+1:N) = 0
227: *> \endverbatim
228: *>
229: *> \param[out] U
230: *> \verbatim
231: *> U is COMPLEX*16 array, dimension (LDU,M)
232: *> If JOBU = 'U', U contains the M-by-M unitary matrix U.
233: *> If JOBU = 'N', U is not referenced.
234: *> \endverbatim
235: *>
236: *> \param[in] LDU
237: *> \verbatim
238: *> LDU is INTEGER
239: *> The leading dimension of the array U. LDU >= max(1,M) if
240: *> JOBU = 'U'; LDU >= 1 otherwise.
241: *> \endverbatim
242: *>
243: *> \param[out] V
244: *> \verbatim
245: *> V is COMPLEX*16 array, dimension (LDV,P)
246: *> If JOBV = 'V', V contains the P-by-P unitary matrix V.
247: *> If JOBV = 'N', V is not referenced.
248: *> \endverbatim
249: *>
250: *> \param[in] LDV
251: *> \verbatim
252: *> LDV is INTEGER
253: *> The leading dimension of the array V. LDV >= max(1,P) if
254: *> JOBV = 'V'; LDV >= 1 otherwise.
255: *> \endverbatim
256: *>
257: *> \param[out] Q
258: *> \verbatim
259: *> Q is COMPLEX*16 array, dimension (LDQ,N)
260: *> If JOBQ = 'Q', Q contains the N-by-N unitary matrix Q.
261: *> If JOBQ = 'N', Q is not referenced.
262: *> \endverbatim
263: *>
264: *> \param[in] LDQ
265: *> \verbatim
266: *> LDQ is INTEGER
267: *> The leading dimension of the array Q. LDQ >= max(1,N) if
268: *> JOBQ = 'Q'; LDQ >= 1 otherwise.
269: *> \endverbatim
270: *>
271: *> \param[out] WORK
272: *> \verbatim
273: *> WORK is COMPLEX*16 array, dimension (MAX(1,LWORK))
274: *> On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
275: *> \endverbatim
276: *>
277: *> \param[in] LWORK
278: *> \verbatim
279: *> LWORK is INTEGER
280: *> The dimension of the array WORK.
281: *>
282: *> If LWORK = -1, then a workspace query is assumed; the routine
283: *> only calculates the optimal size of the WORK array, returns
284: *> this value as the first entry of the WORK array, and no error
285: *> message related to LWORK is issued by XERBLA.
286: *> \endverbatim
287: *>
288: *> \param[out] RWORK
289: *> \verbatim
290: *> RWORK is DOUBLE PRECISION array, dimension (2*N)
291: *> \endverbatim
292: *>
293: *> \param[out] IWORK
294: *> \verbatim
295: *> IWORK is INTEGER array, dimension (N)
296: *> On exit, IWORK stores the sorting information. More
297: *> precisely, the following loop will sort ALPHA
298: *> for I = K+1, min(M,K+L)
299: *> swap ALPHA(I) and ALPHA(IWORK(I))
300: *> endfor
301: *> such that ALPHA(1) >= ALPHA(2) >= ... >= ALPHA(N).
302: *> \endverbatim
303: *>
304: *> \param[out] INFO
305: *> \verbatim
306: *> INFO is INTEGER
307: *> = 0: successful exit.
308: *> < 0: if INFO = -i, the i-th argument had an illegal value.
309: *> > 0: if INFO = 1, the Jacobi-type procedure failed to
310: *> converge. For further details, see subroutine ZTGSJA.
311: *> \endverbatim
312: *
313: *> \par Internal Parameters:
314: * =========================
315: *>
316: *> \verbatim
317: *> TOLA DOUBLE PRECISION
318: *> TOLB DOUBLE PRECISION
319: *> TOLA and TOLB are the thresholds to determine the effective
320: *> rank of (A**H,B**H)**H. Generally, they are set to
321: *> TOLA = MAX(M,N)*norm(A)*MACHEPS,
322: *> TOLB = MAX(P,N)*norm(B)*MACHEPS.
323: *> The size of TOLA and TOLB may affect the size of backward
324: *> errors of the decomposition.
325: *> \endverbatim
326: *
327: * Authors:
328: * ========
329: *
330: *> \author Univ. of Tennessee
331: *> \author Univ. of California Berkeley
332: *> \author Univ. of Colorado Denver
333: *> \author NAG Ltd.
334: *
335: *> \date August 2015
336: *
337: *> \ingroup complex16GEsing
338: *
339: *> \par Contributors:
340: * ==================
341: *>
342: *> Ming Gu and Huan Ren, Computer Science Division, University of
343: *> California at Berkeley, USA
344: *>
345: *
346: *> \par Further Details:
347: * =====================
348: *>
349: *> ZGGSVD3 replaces the deprecated subroutine ZGGSVD.
350: *>
351: * =====================================================================
352: SUBROUTINE ZGGSVD3( JOBU, JOBV, JOBQ, M, N, P, K, L, A, LDA, B,
353: $ LDB, ALPHA, BETA, U, LDU, V, LDV, Q, LDQ,
354: $ WORK, LWORK, RWORK, IWORK, INFO )
355: *
356: * -- LAPACK driver routine (version 3.7.0) --
357: * -- LAPACK is a software package provided by Univ. of Tennessee, --
358: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
359: * August 2015
360: *
361: * .. Scalar Arguments ..
362: CHARACTER JOBQ, JOBU, JOBV
363: INTEGER INFO, K, L, LDA, LDB, LDQ, LDU, LDV, M, N, P,
364: $ LWORK
365: * ..
366: * .. Array Arguments ..
367: INTEGER IWORK( * )
368: DOUBLE PRECISION ALPHA( * ), BETA( * ), RWORK( * )
369: COMPLEX*16 A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
370: $ U( LDU, * ), V( LDV, * ), WORK( * )
371: * ..
372: *
373: * =====================================================================
374: *
375: * .. Local Scalars ..
376: LOGICAL WANTQ, WANTU, WANTV, LQUERY
377: INTEGER I, IBND, ISUB, J, NCYCLE, LWKOPT
378: DOUBLE PRECISION ANORM, BNORM, SMAX, TEMP, TOLA, TOLB, ULP, UNFL
379: * ..
380: * .. External Functions ..
381: LOGICAL LSAME
382: DOUBLE PRECISION DLAMCH, ZLANGE
383: EXTERNAL LSAME, DLAMCH, ZLANGE
384: * ..
385: * .. External Subroutines ..
386: EXTERNAL DCOPY, XERBLA, ZGGSVP3, ZTGSJA
387: * ..
388: * .. Intrinsic Functions ..
389: INTRINSIC MAX, MIN
390: * ..
391: * .. Executable Statements ..
392: *
393: * Decode and test the input parameters
394: *
395: WANTU = LSAME( JOBU, 'U' )
396: WANTV = LSAME( JOBV, 'V' )
397: WANTQ = LSAME( JOBQ, 'Q' )
398: LQUERY = ( LWORK.EQ.-1 )
399: LWKOPT = 1
400: *
401: * Test the input arguments
402: *
403: INFO = 0
404: IF( .NOT.( WANTU .OR. LSAME( JOBU, 'N' ) ) ) THEN
405: INFO = -1
406: ELSE IF( .NOT.( WANTV .OR. LSAME( JOBV, 'N' ) ) ) THEN
407: INFO = -2
408: ELSE IF( .NOT.( WANTQ .OR. LSAME( JOBQ, 'N' ) ) ) THEN
409: INFO = -3
410: ELSE IF( M.LT.0 ) THEN
411: INFO = -4
412: ELSE IF( N.LT.0 ) THEN
413: INFO = -5
414: ELSE IF( P.LT.0 ) THEN
415: INFO = -6
416: ELSE IF( LDA.LT.MAX( 1, M ) ) THEN
417: INFO = -10
418: ELSE IF( LDB.LT.MAX( 1, P ) ) THEN
419: INFO = -12
420: ELSE IF( LDU.LT.1 .OR. ( WANTU .AND. LDU.LT.M ) ) THEN
421: INFO = -16
422: ELSE IF( LDV.LT.1 .OR. ( WANTV .AND. LDV.LT.P ) ) THEN
423: INFO = -18
424: ELSE IF( LDQ.LT.1 .OR. ( WANTQ .AND. LDQ.LT.N ) ) THEN
425: INFO = -20
426: ELSE IF( LWORK.LT.1 .AND. .NOT.LQUERY ) THEN
427: INFO = -24
428: END IF
429: *
430: * Compute workspace
431: *
432: IF( INFO.EQ.0 ) THEN
433: CALL ZGGSVP3( JOBU, JOBV, JOBQ, M, P, N, A, LDA, B, LDB, TOLA,
434: $ TOLB, K, L, U, LDU, V, LDV, Q, LDQ, IWORK, RWORK,
435: $ WORK, WORK, -1, INFO )
436: LWKOPT = N + INT( WORK( 1 ) )
437: LWKOPT = MAX( 2*N, LWKOPT )
438: LWKOPT = MAX( 1, LWKOPT )
439: WORK( 1 ) = DCMPLX( LWKOPT )
440: END IF
441: *
442: IF( INFO.NE.0 ) THEN
443: CALL XERBLA( 'ZGGSVD3', -INFO )
444: RETURN
445: END IF
446: IF( LQUERY ) THEN
447: RETURN
448: ENDIF
449: *
450: * Compute the Frobenius norm of matrices A and B
451: *
452: ANORM = ZLANGE( '1', M, N, A, LDA, RWORK )
453: BNORM = ZLANGE( '1', P, N, B, LDB, RWORK )
454: *
455: * Get machine precision and set up threshold for determining
456: * the effective numerical rank of the matrices A and B.
457: *
458: ULP = DLAMCH( 'Precision' )
459: UNFL = DLAMCH( 'Safe Minimum' )
460: TOLA = MAX( M, N )*MAX( ANORM, UNFL )*ULP
461: TOLB = MAX( P, N )*MAX( BNORM, UNFL )*ULP
462: *
463: CALL ZGGSVP3( JOBU, JOBV, JOBQ, M, P, N, A, LDA, B, LDB, TOLA,
464: $ TOLB, K, L, U, LDU, V, LDV, Q, LDQ, IWORK, RWORK,
465: $ WORK, WORK( N+1 ), LWORK-N, INFO )
466: *
467: * Compute the GSVD of two upper "triangular" matrices
468: *
469: CALL ZTGSJA( JOBU, JOBV, JOBQ, M, P, N, K, L, A, LDA, B, LDB,
470: $ TOLA, TOLB, ALPHA, BETA, U, LDU, V, LDV, Q, LDQ,
471: $ WORK, NCYCLE, INFO )
472: *
473: * Sort the singular values and store the pivot indices in IWORK
474: * Copy ALPHA to RWORK, then sort ALPHA in RWORK
475: *
476: CALL DCOPY( N, ALPHA, 1, RWORK, 1 )
477: IBND = MIN( L, M-K )
478: DO 20 I = 1, IBND
479: *
480: * Scan for largest ALPHA(K+I)
481: *
482: ISUB = I
483: SMAX = RWORK( K+I )
484: DO 10 J = I + 1, IBND
485: TEMP = RWORK( K+J )
486: IF( TEMP.GT.SMAX ) THEN
487: ISUB = J
488: SMAX = TEMP
489: END IF
490: 10 CONTINUE
491: IF( ISUB.NE.I ) THEN
492: RWORK( K+ISUB ) = RWORK( K+I )
493: RWORK( K+I ) = SMAX
494: IWORK( K+I ) = K + ISUB
495: ELSE
496: IWORK( K+I ) = K + I
497: END IF
498: 20 CONTINUE
499: *
500: WORK( 1 ) = DCMPLX( LWKOPT )
501: RETURN
502: *
503: * End of ZGGSVD3
504: *
505: END
CVSweb interface <joel.bertrand@systella.fr>