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