LAPACK 3.3.0
|
00001 SUBROUTINE ZGGSVD( JOBU, JOBV, JOBQ, M, N, P, K, L, A, LDA, B, 00002 $ LDB, ALPHA, BETA, U, LDU, V, LDV, Q, LDQ, WORK, 00003 $ RWORK, IWORK, INFO ) 00004 * 00005 * -- LAPACK driver routine (version 3.2) -- 00006 * -- LAPACK is a software package provided by Univ. of Tennessee, -- 00007 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 00008 * November 2006 00009 * 00010 * .. Scalar Arguments .. 00011 CHARACTER JOBQ, JOBU, JOBV 00012 INTEGER INFO, K, L, LDA, LDB, LDQ, LDU, LDV, M, N, P 00013 * .. 00014 * .. Array Arguments .. 00015 INTEGER IWORK( * ) 00016 DOUBLE PRECISION ALPHA( * ), BETA( * ), RWORK( * ) 00017 COMPLEX*16 A( LDA, * ), B( LDB, * ), Q( LDQ, * ), 00018 $ U( LDU, * ), V( LDV, * ), WORK( * ) 00019 * .. 00020 * 00021 * Purpose 00022 * ======= 00023 * 00024 * ZGGSVD computes the generalized singular value decomposition (GSVD) 00025 * of an M-by-N complex matrix A and P-by-N complex matrix B: 00026 * 00027 * U'*A*Q = D1*( 0 R ), V'*B*Q = D2*( 0 R ) 00028 * 00029 * where U, V and Q are unitary matrices, and Z' means the conjugate 00030 * transpose of Z. Let K+L = the effective numerical rank of the 00031 * matrix (A',B')', then R is a (K+L)-by-(K+L) nonsingular upper 00032 * triangular matrix, D1 and D2 are M-by-(K+L) and P-by-(K+L) "diagonal" 00033 * matrices and of the following structures, respectively: 00034 * 00035 * If M-K-L >= 0, 00036 * 00037 * K L 00038 * D1 = K ( I 0 ) 00039 * L ( 0 C ) 00040 * M-K-L ( 0 0 ) 00041 * 00042 * K L 00043 * D2 = L ( 0 S ) 00044 * P-L ( 0 0 ) 00045 * 00046 * N-K-L K L 00047 * ( 0 R ) = K ( 0 R11 R12 ) 00048 * L ( 0 0 R22 ) 00049 * where 00050 * 00051 * C = diag( ALPHA(K+1), ... , ALPHA(K+L) ), 00052 * S = diag( BETA(K+1), ... , BETA(K+L) ), 00053 * C**2 + S**2 = I. 00054 * 00055 * R is stored in A(1:K+L,N-K-L+1:N) on exit. 00056 * 00057 * If M-K-L < 0, 00058 * 00059 * K M-K K+L-M 00060 * D1 = K ( I 0 0 ) 00061 * M-K ( 0 C 0 ) 00062 * 00063 * K M-K K+L-M 00064 * D2 = M-K ( 0 S 0 ) 00065 * K+L-M ( 0 0 I ) 00066 * P-L ( 0 0 0 ) 00067 * 00068 * N-K-L K M-K K+L-M 00069 * ( 0 R ) = K ( 0 R11 R12 R13 ) 00070 * M-K ( 0 0 R22 R23 ) 00071 * K+L-M ( 0 0 0 R33 ) 00072 * 00073 * where 00074 * 00075 * C = diag( ALPHA(K+1), ... , ALPHA(M) ), 00076 * S = diag( BETA(K+1), ... , BETA(M) ), 00077 * C**2 + S**2 = I. 00078 * 00079 * (R11 R12 R13 ) is stored in A(1:M, N-K-L+1:N), and R33 is stored 00080 * ( 0 R22 R23 ) 00081 * in B(M-K+1:L,N+M-K-L+1:N) on exit. 00082 * 00083 * The routine computes C, S, R, and optionally the unitary 00084 * transformation matrices U, V and Q. 00085 * 00086 * In particular, if B is an N-by-N nonsingular matrix, then the GSVD of 00087 * A and B implicitly gives the SVD of A*inv(B): 00088 * A*inv(B) = U*(D1*inv(D2))*V'. 00089 * If ( A',B')' has orthnormal columns, then the GSVD of A and B is also 00090 * equal to the CS decomposition of A and B. Furthermore, the GSVD can 00091 * be used to derive the solution of the eigenvalue problem: 00092 * A'*A x = lambda* B'*B x. 00093 * In some literature, the GSVD of A and B is presented in the form 00094 * U'*A*X = ( 0 D1 ), V'*B*X = ( 0 D2 ) 00095 * where U and V are orthogonal and X is nonsingular, and D1 and D2 are 00096 * ``diagonal''. The former GSVD form can be converted to the latter 00097 * form by taking the nonsingular matrix X as 00098 * 00099 * X = Q*( I 0 ) 00100 * ( 0 inv(R) ) 00101 * 00102 * Arguments 00103 * ========= 00104 * 00105 * JOBU (input) CHARACTER*1 00106 * = 'U': Unitary matrix U is computed; 00107 * = 'N': U is not computed. 00108 * 00109 * JOBV (input) CHARACTER*1 00110 * = 'V': Unitary matrix V is computed; 00111 * = 'N': V is not computed. 00112 * 00113 * JOBQ (input) CHARACTER*1 00114 * = 'Q': Unitary matrix Q is computed; 00115 * = 'N': Q is not computed. 00116 * 00117 * M (input) INTEGER 00118 * The number of rows of the matrix A. M >= 0. 00119 * 00120 * N (input) INTEGER 00121 * The number of columns of the matrices A and B. N >= 0. 00122 * 00123 * P (input) INTEGER 00124 * The number of rows of the matrix B. P >= 0. 00125 * 00126 * K (output) INTEGER 00127 * L (output) INTEGER 00128 * On exit, K and L specify the dimension of the subblocks 00129 * described in Purpose. 00130 * K + L = effective numerical rank of (A',B')'. 00131 * 00132 * A (input/output) COMPLEX*16 array, dimension (LDA,N) 00133 * On entry, the M-by-N matrix A. 00134 * On exit, A contains the triangular matrix R, or part of R. 00135 * See Purpose for details. 00136 * 00137 * LDA (input) INTEGER 00138 * The leading dimension of the array A. LDA >= max(1,M). 00139 * 00140 * B (input/output) COMPLEX*16 array, dimension (LDB,N) 00141 * On entry, the P-by-N matrix B. 00142 * On exit, B contains part of the triangular matrix R if 00143 * M-K-L < 0. See Purpose for details. 00144 * 00145 * LDB (input) INTEGER 00146 * The leading dimension of the array B. LDB >= max(1,P). 00147 * 00148 * ALPHA (output) DOUBLE PRECISION array, dimension (N) 00149 * BETA (output) DOUBLE PRECISION array, dimension (N) 00150 * On exit, ALPHA and BETA contain the generalized singular 00151 * value pairs of A and B; 00152 * ALPHA(1:K) = 1, 00153 * BETA(1:K) = 0, 00154 * and if M-K-L >= 0, 00155 * ALPHA(K+1:K+L) = C, 00156 * BETA(K+1:K+L) = S, 00157 * or if M-K-L < 0, 00158 * ALPHA(K+1:M)= C, ALPHA(M+1:K+L)= 0 00159 * BETA(K+1:M) = S, BETA(M+1:K+L) = 1 00160 * and 00161 * ALPHA(K+L+1:N) = 0 00162 * BETA(K+L+1:N) = 0 00163 * 00164 * U (output) COMPLEX*16 array, dimension (LDU,M) 00165 * If JOBU = 'U', U contains the M-by-M unitary matrix U. 00166 * If JOBU = 'N', U is not referenced. 00167 * 00168 * LDU (input) INTEGER 00169 * The leading dimension of the array U. LDU >= max(1,M) if 00170 * JOBU = 'U'; LDU >= 1 otherwise. 00171 * 00172 * V (output) COMPLEX*16 array, dimension (LDV,P) 00173 * If JOBV = 'V', V contains the P-by-P unitary matrix V. 00174 * If JOBV = 'N', V is not referenced. 00175 * 00176 * LDV (input) INTEGER 00177 * The leading dimension of the array V. LDV >= max(1,P) if 00178 * JOBV = 'V'; LDV >= 1 otherwise. 00179 * 00180 * Q (output) COMPLEX*16 array, dimension (LDQ,N) 00181 * If JOBQ = 'Q', Q contains the N-by-N unitary matrix Q. 00182 * If JOBQ = 'N', Q is not referenced. 00183 * 00184 * LDQ (input) INTEGER 00185 * The leading dimension of the array Q. LDQ >= max(1,N) if 00186 * JOBQ = 'Q'; LDQ >= 1 otherwise. 00187 * 00188 * WORK (workspace) COMPLEX*16 array, dimension (max(3*N,M,P)+N) 00189 * 00190 * RWORK (workspace) DOUBLE PRECISION array, dimension (2*N) 00191 * 00192 * IWORK (workspace/output) INTEGER array, dimension (N) 00193 * On exit, IWORK stores the sorting information. More 00194 * precisely, the following loop will sort ALPHA 00195 * for I = K+1, min(M,K+L) 00196 * swap ALPHA(I) and ALPHA(IWORK(I)) 00197 * endfor 00198 * such that ALPHA(1) >= ALPHA(2) >= ... >= ALPHA(N). 00199 * 00200 * INFO (output) INTEGER 00201 * = 0: successful exit. 00202 * < 0: if INFO = -i, the i-th argument had an illegal value. 00203 * > 0: if INFO = 1, the Jacobi-type procedure failed to 00204 * converge. For further details, see subroutine ZTGSJA. 00205 * 00206 * Internal Parameters 00207 * =================== 00208 * 00209 * TOLA DOUBLE PRECISION 00210 * TOLB DOUBLE PRECISION 00211 * TOLA and TOLB are the thresholds to determine the effective 00212 * rank of (A',B')'. Generally, they are set to 00213 * TOLA = MAX(M,N)*norm(A)*MAZHEPS, 00214 * TOLB = MAX(P,N)*norm(B)*MAZHEPS. 00215 * The size of TOLA and TOLB may affect the size of backward 00216 * errors of the decomposition. 00217 * 00218 * Further Details 00219 * =============== 00220 * 00221 * 2-96 Based on modifications by 00222 * Ming Gu and Huan Ren, Computer Science Division, University of 00223 * California at Berkeley, USA 00224 * 00225 * ===================================================================== 00226 * 00227 * .. Local Scalars .. 00228 LOGICAL WANTQ, WANTU, WANTV 00229 INTEGER I, IBND, ISUB, J, NCYCLE 00230 DOUBLE PRECISION ANORM, BNORM, SMAX, TEMP, TOLA, TOLB, ULP, UNFL 00231 * .. 00232 * .. External Functions .. 00233 LOGICAL LSAME 00234 DOUBLE PRECISION DLAMCH, ZLANGE 00235 EXTERNAL LSAME, DLAMCH, ZLANGE 00236 * .. 00237 * .. External Subroutines .. 00238 EXTERNAL DCOPY, XERBLA, ZGGSVP, ZTGSJA 00239 * .. 00240 * .. Intrinsic Functions .. 00241 INTRINSIC MAX, MIN 00242 * .. 00243 * .. Executable Statements .. 00244 * 00245 * Decode and test the input parameters 00246 * 00247 WANTU = LSAME( JOBU, 'U' ) 00248 WANTV = LSAME( JOBV, 'V' ) 00249 WANTQ = LSAME( JOBQ, 'Q' ) 00250 * 00251 INFO = 0 00252 IF( .NOT.( WANTU .OR. LSAME( JOBU, 'N' ) ) ) THEN 00253 INFO = -1 00254 ELSE IF( .NOT.( WANTV .OR. LSAME( JOBV, 'N' ) ) ) THEN 00255 INFO = -2 00256 ELSE IF( .NOT.( WANTQ .OR. LSAME( JOBQ, 'N' ) ) ) THEN 00257 INFO = -3 00258 ELSE IF( M.LT.0 ) THEN 00259 INFO = -4 00260 ELSE IF( N.LT.0 ) THEN 00261 INFO = -5 00262 ELSE IF( P.LT.0 ) THEN 00263 INFO = -6 00264 ELSE IF( LDA.LT.MAX( 1, M ) ) THEN 00265 INFO = -10 00266 ELSE IF( LDB.LT.MAX( 1, P ) ) THEN 00267 INFO = -12 00268 ELSE IF( LDU.LT.1 .OR. ( WANTU .AND. LDU.LT.M ) ) THEN 00269 INFO = -16 00270 ELSE IF( LDV.LT.1 .OR. ( WANTV .AND. LDV.LT.P ) ) THEN 00271 INFO = -18 00272 ELSE IF( LDQ.LT.1 .OR. ( WANTQ .AND. LDQ.LT.N ) ) THEN 00273 INFO = -20 00274 END IF 00275 IF( INFO.NE.0 ) THEN 00276 CALL XERBLA( 'ZGGSVD', -INFO ) 00277 RETURN 00278 END IF 00279 * 00280 * Compute the Frobenius norm of matrices A and B 00281 * 00282 ANORM = ZLANGE( '1', M, N, A, LDA, RWORK ) 00283 BNORM = ZLANGE( '1', P, N, B, LDB, RWORK ) 00284 * 00285 * Get machine precision and set up threshold for determining 00286 * the effective numerical rank of the matrices A and B. 00287 * 00288 ULP = DLAMCH( 'Precision' ) 00289 UNFL = DLAMCH( 'Safe Minimum' ) 00290 TOLA = MAX( M, N )*MAX( ANORM, UNFL )*ULP 00291 TOLB = MAX( P, N )*MAX( BNORM, UNFL )*ULP 00292 * 00293 CALL ZGGSVP( JOBU, JOBV, JOBQ, M, P, N, A, LDA, B, LDB, TOLA, 00294 $ TOLB, K, L, U, LDU, V, LDV, Q, LDQ, IWORK, RWORK, 00295 $ WORK, WORK( N+1 ), INFO ) 00296 * 00297 * Compute the GSVD of two upper "triangular" matrices 00298 * 00299 CALL ZTGSJA( JOBU, JOBV, JOBQ, M, P, N, K, L, A, LDA, B, LDB, 00300 $ TOLA, TOLB, ALPHA, BETA, U, LDU, V, LDV, Q, LDQ, 00301 $ WORK, NCYCLE, INFO ) 00302 * 00303 * Sort the singular values and store the pivot indices in IWORK 00304 * Copy ALPHA to RWORK, then sort ALPHA in RWORK 00305 * 00306 CALL DCOPY( N, ALPHA, 1, RWORK, 1 ) 00307 IBND = MIN( L, M-K ) 00308 DO 20 I = 1, IBND 00309 * 00310 * Scan for largest ALPHA(K+I) 00311 * 00312 ISUB = I 00313 SMAX = RWORK( K+I ) 00314 DO 10 J = I + 1, IBND 00315 TEMP = RWORK( K+J ) 00316 IF( TEMP.GT.SMAX ) THEN 00317 ISUB = J 00318 SMAX = TEMP 00319 END IF 00320 10 CONTINUE 00321 IF( ISUB.NE.I ) THEN 00322 RWORK( K+ISUB ) = RWORK( K+I ) 00323 RWORK( K+I ) = SMAX 00324 IWORK( K+I ) = K + ISUB 00325 ELSE 00326 IWORK( K+I ) = K + I 00327 END IF 00328 20 CONTINUE 00329 * 00330 RETURN 00331 * 00332 * End of ZGGSVD 00333 * 00334 END