LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sggsvp.f
Go to the documentation of this file.
1*> \brief \b SGGSVP
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download SGGSVP + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sggsvp.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sggsvp.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sggsvp.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE SGGSVP( JOBU, JOBV, JOBQ, M, P, N, A, LDA, B, LDB,
22* TOLA, TOLB, K, L, U, LDU, V, LDV, Q, LDQ,
23* IWORK, TAU, WORK, INFO )
24*
25* .. Scalar Arguments ..
26* CHARACTER JOBQ, JOBU, JOBV
27* INTEGER INFO, K, L, LDA, LDB, LDQ, LDU, LDV, M, N, P
28* REAL TOLA, TOLB
29* ..
30* .. Array Arguments ..
31* INTEGER IWORK( * )
32* REAL A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
33* $ TAU( * ), U( LDU, * ), V( LDV, * ), WORK( * )
34* ..
35*
36*
37*> \par Purpose:
38* =============
39*>
40*> \verbatim
41*>
42*> This routine is deprecated and has been replaced by routine SGGSVP3.
43*>
44*> SGGSVP computes orthogonal matrices U, V and Q such that
45*>
46*> N-K-L K L
47*> U**T*A*Q = K ( 0 A12 A13 ) if M-K-L >= 0;
48*> L ( 0 0 A23 )
49*> M-K-L ( 0 0 0 )
50*>
51*> N-K-L K L
52*> = K ( 0 A12 A13 ) if M-K-L < 0;
53*> M-K ( 0 0 A23 )
54*>
55*> N-K-L K L
56*> V**T*B*Q = L ( 0 0 B13 )
57*> P-L ( 0 0 0 )
58*>
59*> where the K-by-K matrix A12 and L-by-L matrix B13 are nonsingular
60*> upper triangular; A23 is L-by-L upper triangular if M-K-L >= 0,
61*> otherwise A23 is (M-K)-by-L upper trapezoidal. K+L = the effective
62*> numerical rank of the (M+P)-by-N matrix (A**T,B**T)**T.
63*>
64*> This decomposition is the preprocessing step for computing the
65*> Generalized Singular Value Decomposition (GSVD), see subroutine
66*> SGGSVD.
67*> \endverbatim
68*
69* Arguments:
70* ==========
71*
72*> \param[in] JOBU
73*> \verbatim
74*> JOBU is CHARACTER*1
75*> = 'U': Orthogonal matrix U is computed;
76*> = 'N': U is not computed.
77*> \endverbatim
78*>
79*> \param[in] JOBV
80*> \verbatim
81*> JOBV is CHARACTER*1
82*> = 'V': Orthogonal matrix V is computed;
83*> = 'N': V is not computed.
84*> \endverbatim
85*>
86*> \param[in] JOBQ
87*> \verbatim
88*> JOBQ is CHARACTER*1
89*> = 'Q': Orthogonal matrix Q is computed;
90*> = 'N': Q is not computed.
91*> \endverbatim
92*>
93*> \param[in] M
94*> \verbatim
95*> M is INTEGER
96*> The number of rows of the matrix A. M >= 0.
97*> \endverbatim
98*>
99*> \param[in] P
100*> \verbatim
101*> P is INTEGER
102*> The number of rows of the matrix B. P >= 0.
103*> \endverbatim
104*>
105*> \param[in] N
106*> \verbatim
107*> N is INTEGER
108*> The number of columns of the matrices A and B. N >= 0.
109*> \endverbatim
110*>
111*> \param[in,out] A
112*> \verbatim
113*> A is REAL array, dimension (LDA,N)
114*> On entry, the M-by-N matrix A.
115*> On exit, A contains the triangular (or trapezoidal) matrix
116*> described in the Purpose section.
117*> \endverbatim
118*>
119*> \param[in] LDA
120*> \verbatim
121*> LDA is INTEGER
122*> The leading dimension of the array A. LDA >= max(1,M).
123*> \endverbatim
124*>
125*> \param[in,out] B
126*> \verbatim
127*> B is REAL array, dimension (LDB,N)
128*> On entry, the P-by-N matrix B.
129*> On exit, B contains the triangular matrix described in
130*> the Purpose section.
131*> \endverbatim
132*>
133*> \param[in] LDB
134*> \verbatim
135*> LDB is INTEGER
136*> The leading dimension of the array B. LDB >= max(1,P).
137*> \endverbatim
138*>
139*> \param[in] TOLA
140*> \verbatim
141*> TOLA is REAL
142*> \endverbatim
143*>
144*> \param[in] TOLB
145*> \verbatim
146*> TOLB is REAL
147*>
148*> TOLA and TOLB are the thresholds to determine the effective
149*> numerical rank of matrix B and a subblock of A. Generally,
150*> they are set to
151*> TOLA = MAX(M,N)*norm(A)*MACHEPS,
152*> TOLB = MAX(P,N)*norm(B)*MACHEPS.
153*> The size of TOLA and TOLB may affect the size of backward
154*> errors of the decomposition.
155*> \endverbatim
156*>
157*> \param[out] K
158*> \verbatim
159*> K is INTEGER
160*> \endverbatim
161*>
162*> \param[out] L
163*> \verbatim
164*> L is INTEGER
165*>
166*> On exit, K and L specify the dimension of the subblocks
167*> described in Purpose section.
168*> K + L = effective numerical rank of (A**T,B**T)**T.
169*> \endverbatim
170*>
171*> \param[out] U
172*> \verbatim
173*> U is REAL array, dimension (LDU,M)
174*> If JOBU = 'U', U contains the orthogonal matrix U.
175*> If JOBU = 'N', U is not referenced.
176*> \endverbatim
177*>
178*> \param[in] LDU
179*> \verbatim
180*> LDU is INTEGER
181*> The leading dimension of the array U. LDU >= max(1,M) if
182*> JOBU = 'U'; LDU >= 1 otherwise.
183*> \endverbatim
184*>
185*> \param[out] V
186*> \verbatim
187*> V is REAL array, dimension (LDV,P)
188*> If JOBV = 'V', V contains the orthogonal matrix V.
189*> If JOBV = 'N', V is not referenced.
190*> \endverbatim
191*>
192*> \param[in] LDV
193*> \verbatim
194*> LDV is INTEGER
195*> The leading dimension of the array V. LDV >= max(1,P) if
196*> JOBV = 'V'; LDV >= 1 otherwise.
197*> \endverbatim
198*>
199*> \param[out] Q
200*> \verbatim
201*> Q is REAL array, dimension (LDQ,N)
202*> If JOBQ = 'Q', Q contains the orthogonal matrix Q.
203*> If JOBQ = 'N', Q is not referenced.
204*> \endverbatim
205*>
206*> \param[in] LDQ
207*> \verbatim
208*> LDQ is INTEGER
209*> The leading dimension of the array Q. LDQ >= max(1,N) if
210*> JOBQ = 'Q'; LDQ >= 1 otherwise.
211*> \endverbatim
212*>
213*> \param[out] IWORK
214*> \verbatim
215*> IWORK is INTEGER array, dimension (N)
216*> \endverbatim
217*>
218*> \param[out] TAU
219*> \verbatim
220*> TAU is REAL array, dimension (N)
221*> \endverbatim
222*>
223*> \param[out] WORK
224*> \verbatim
225*> WORK is REAL array, dimension (max(3*N,M,P))
226*> \endverbatim
227*>
228*> \param[out] INFO
229*> \verbatim
230*> INFO is INTEGER
231*> = 0: successful exit
232*> < 0: if INFO = -i, the i-th argument had an illegal value.
233*> \endverbatim
234*
235* Authors:
236* ========
237*
238*> \author Univ. of Tennessee
239*> \author Univ. of California Berkeley
240*> \author Univ. of Colorado Denver
241*> \author NAG Ltd.
242*
243*> \ingroup realOTHERcomputational
244*
245*> \par Further Details:
246* =====================
247*>
248*> The subroutine uses LAPACK subroutine SGEQPF for the QR factorization
249*> with column pivoting to detect the effective numerical rank of the
250*> a matrix. It may be replaced by a better rank determination strategy.
251*>
252* =====================================================================
253 SUBROUTINE sggsvp( JOBU, JOBV, JOBQ, M, P, N, A, LDA, B, LDB,
254 $ TOLA, TOLB, K, L, U, LDU, V, LDV, Q, LDQ,
255 $ IWORK, TAU, WORK, INFO )
256*
257* -- LAPACK computational routine --
258* -- LAPACK is a software package provided by Univ. of Tennessee, --
259* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
260*
261* .. Scalar Arguments ..
262 CHARACTER JOBQ, JOBU, JOBV
263 INTEGER INFO, K, L, LDA, LDB, LDQ, LDU, LDV, M, N, P
264 REAL TOLA, TOLB
265* ..
266* .. Array Arguments ..
267 INTEGER IWORK( * )
268 REAL A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
269 $ tau( * ), u( ldu, * ), v( ldv, * ), work( * )
270* ..
271*
272* =====================================================================
273*
274* .. Parameters ..
275 REAL ZERO, ONE
276 PARAMETER ( ZERO = 0.0e+0, one = 1.0e+0 )
277* ..
278* .. Local Scalars ..
279 LOGICAL FORWRD, WANTQ, WANTU, WANTV
280 INTEGER I, J
281* ..
282* .. External Functions ..
283 LOGICAL LSAME
284 EXTERNAL LSAME
285* ..
286* .. External Subroutines ..
287 EXTERNAL sgeqpf, sgeqr2, sgerq2, slacpy, slapmt, slaset,
289* ..
290* .. Intrinsic Functions ..
291 INTRINSIC abs, max, min
292* ..
293* .. Executable Statements ..
294*
295* Test the input parameters
296*
297 wantu = lsame( jobu, 'U' )
298 wantv = lsame( jobv, 'V' )
299 wantq = lsame( jobq, 'Q' )
300 forwrd = .true.
301*
302 info = 0
303 IF( .NOT.( wantu .OR. lsame( jobu, 'N' ) ) ) THEN
304 info = -1
305 ELSE IF( .NOT.( wantv .OR. lsame( jobv, 'N' ) ) ) THEN
306 info = -2
307 ELSE IF( .NOT.( wantq .OR. lsame( jobq, 'N' ) ) ) THEN
308 info = -3
309 ELSE IF( m.LT.0 ) THEN
310 info = -4
311 ELSE IF( p.LT.0 ) THEN
312 info = -5
313 ELSE IF( n.LT.0 ) THEN
314 info = -6
315 ELSE IF( lda.LT.max( 1, m ) ) THEN
316 info = -8
317 ELSE IF( ldb.LT.max( 1, p ) ) THEN
318 info = -10
319 ELSE IF( ldu.LT.1 .OR. ( wantu .AND. ldu.LT.m ) ) THEN
320 info = -16
321 ELSE IF( ldv.LT.1 .OR. ( wantv .AND. ldv.LT.p ) ) THEN
322 info = -18
323 ELSE IF( ldq.LT.1 .OR. ( wantq .AND. ldq.LT.n ) ) THEN
324 info = -20
325 END IF
326 IF( info.NE.0 ) THEN
327 CALL xerbla( 'SGGSVP', -info )
328 RETURN
329 END IF
330*
331* QR with column pivoting of B: B*P = V*( S11 S12 )
332* ( 0 0 )
333*
334 DO 10 i = 1, n
335 iwork( i ) = 0
336 10 CONTINUE
337 CALL sgeqpf( p, n, b, ldb, iwork, tau, work, info )
338*
339* Update A := A*P
340*
341 CALL slapmt( forwrd, m, n, a, lda, iwork )
342*
343* Determine the effective rank of matrix B.
344*
345 l = 0
346 DO 20 i = 1, min( p, n )
347 IF( abs( b( i, i ) ).GT.tolb )
348 $ l = l + 1
349 20 CONTINUE
350*
351 IF( wantv ) THEN
352*
353* Copy the details of V, and form V.
354*
355 CALL slaset( 'Full', p, p, zero, zero, v, ldv )
356 IF( p.GT.1 )
357 $ CALL slacpy( 'Lower', p-1, n, b( 2, 1 ), ldb, v( 2, 1 ),
358 $ ldv )
359 CALL sorg2r( p, p, min( p, n ), v, ldv, tau, work, info )
360 END IF
361*
362* Clean up B
363*
364 DO 40 j = 1, l - 1
365 DO 30 i = j + 1, l
366 b( i, j ) = zero
367 30 CONTINUE
368 40 CONTINUE
369 IF( p.GT.l )
370 $ CALL slaset( 'Full', p-l, n, zero, zero, b( l+1, 1 ), ldb )
371*
372 IF( wantq ) THEN
373*
374* Set Q = I and Update Q := Q*P
375*
376 CALL slaset( 'Full', n, n, zero, one, q, ldq )
377 CALL slapmt( forwrd, n, n, q, ldq, iwork )
378 END IF
379*
380 IF( p.GE.l .AND. n.NE.l ) THEN
381*
382* RQ factorization of (S11 S12): ( S11 S12 ) = ( 0 S12 )*Z
383*
384 CALL sgerq2( l, n, b, ldb, tau, work, info )
385*
386* Update A := A*Z**T
387*
388 CALL sormr2( 'Right', 'Transpose', m, n, l, b, ldb, tau, a,
389 $ lda, work, info )
390*
391 IF( wantq ) THEN
392*
393* Update Q := Q*Z**T
394*
395 CALL sormr2( 'Right', 'Transpose', n, n, l, b, ldb, tau, q,
396 $ ldq, work, info )
397 END IF
398*
399* Clean up B
400*
401 CALL slaset( 'Full', l, n-l, zero, zero, b, ldb )
402 DO 60 j = n - l + 1, n
403 DO 50 i = j - n + l + 1, l
404 b( i, j ) = zero
405 50 CONTINUE
406 60 CONTINUE
407*
408 END IF
409*
410* Let N-L L
411* A = ( A11 A12 ) M,
412*
413* then the following does the complete QR decomposition of A11:
414*
415* A11 = U*( 0 T12 )*P1**T
416* ( 0 0 )
417*
418 DO 70 i = 1, n - l
419 iwork( i ) = 0
420 70 CONTINUE
421 CALL sgeqpf( m, n-l, a, lda, iwork, tau, work, info )
422*
423* Determine the effective rank of A11
424*
425 k = 0
426 DO 80 i = 1, min( m, n-l )
427 IF( abs( a( i, i ) ).GT.tola )
428 $ k = k + 1
429 80 CONTINUE
430*
431* Update A12 := U**T*A12, where A12 = A( 1:M, N-L+1:N )
432*
433 CALL sorm2r( 'Left', 'Transpose', m, l, min( m, n-l ), a, lda,
434 $ tau, a( 1, n-l+1 ), lda, work, info )
435*
436 IF( wantu ) THEN
437*
438* Copy the details of U, and form U
439*
440 CALL slaset( 'Full', m, m, zero, zero, u, ldu )
441 IF( m.GT.1 )
442 $ CALL slacpy( 'Lower', m-1, n-l, a( 2, 1 ), lda, u( 2, 1 ),
443 $ ldu )
444 CALL sorg2r( m, m, min( m, n-l ), u, ldu, tau, work, info )
445 END IF
446*
447 IF( wantq ) THEN
448*
449* Update Q( 1:N, 1:N-L ) = Q( 1:N, 1:N-L )*P1
450*
451 CALL slapmt( forwrd, n, n-l, q, ldq, iwork )
452 END IF
453*
454* Clean up A: set the strictly lower triangular part of
455* A(1:K, 1:K) = 0, and A( K+1:M, 1:N-L ) = 0.
456*
457 DO 100 j = 1, k - 1
458 DO 90 i = j + 1, k
459 a( i, j ) = zero
460 90 CONTINUE
461 100 CONTINUE
462 IF( m.GT.k )
463 $ CALL slaset( 'Full', m-k, n-l, zero, zero, a( k+1, 1 ), lda )
464*
465 IF( n-l.GT.k ) THEN
466*
467* RQ factorization of ( T11 T12 ) = ( 0 T12 )*Z1
468*
469 CALL sgerq2( k, n-l, a, lda, tau, work, info )
470*
471 IF( wantq ) THEN
472*
473* Update Q( 1:N,1:N-L ) = Q( 1:N,1:N-L )*Z1**T
474*
475 CALL sormr2( 'Right', 'Transpose', n, n-l, k, a, lda, tau,
476 $ q, ldq, work, info )
477 END IF
478*
479* Clean up A
480*
481 CALL slaset( 'Full', k, n-l-k, zero, zero, a, lda )
482 DO 120 j = n - l - k + 1, n - l
483 DO 110 i = j - n + l + k + 1, k
484 a( i, j ) = zero
485 110 CONTINUE
486 120 CONTINUE
487*
488 END IF
489*
490 IF( m.GT.k ) THEN
491*
492* QR factorization of A( K+1:M,N-L+1:N )
493*
494 CALL sgeqr2( m-k, l, a( k+1, n-l+1 ), lda, tau, work, info )
495*
496 IF( wantu ) THEN
497*
498* Update U(:,K+1:M) := U(:,K+1:M)*U1
499*
500 CALL sorm2r( 'Right', 'No transpose', m, m-k, min( m-k, l ),
501 $ a( k+1, n-l+1 ), lda, tau, u( 1, k+1 ), ldu,
502 $ work, info )
503 END IF
504*
505* Clean up
506*
507 DO 140 j = n - l + 1, n
508 DO 130 i = j - n + k + l + 1, m
509 a( i, j ) = zero
510 130 CONTINUE
511 140 CONTINUE
512*
513 END IF
514*
515 RETURN
516*
517* End of SGGSVP
518*
519 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine sgeqr2(m, n, a, lda, tau, work, info)
SGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm.
Definition sgeqr2.f:130
subroutine sgerq2(m, n, a, lda, tau, work, info)
SGERQ2 computes the RQ factorization of a general rectangular matrix using an unblocked algorithm.
Definition sgerq2.f:123
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
subroutine slapmt(forwrd, m, n, x, ldx, k)
SLAPMT performs a forward or backward permutation of the columns of a matrix.
Definition slapmt.f:104
subroutine slaset(uplo, m, n, alpha, beta, a, lda)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition slaset.f:110
subroutine sorg2r(m, n, k, a, lda, tau, work, info)
SORG2R generates all or part of the orthogonal matrix Q from a QR factorization determined by sgeqrf ...
Definition sorg2r.f:114
subroutine sorm2r(side, trans, m, n, k, a, lda, tau, c, ldc, work, info)
SORM2R multiplies a general matrix by the orthogonal matrix from a QR factorization determined by sge...
Definition sorm2r.f:159
subroutine sormr2(side, trans, m, n, k, a, lda, tau, c, ldc, work, info)
SORMR2 multiplies a general matrix by the orthogonal matrix from a RQ factorization determined by sge...
Definition sormr2.f:159
subroutine sgeqpf(m, n, a, lda, jpvt, tau, work, info)
SGEQPF
Definition sgeqpf.f:142
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:256