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