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