LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ sggsvp()

subroutine sggsvp ( character jobu,
character jobv,
character jobq,
integer m,
integer p,
integer n,
real, dimension( lda, * ) a,
integer lda,
real, dimension( ldb, * ) b,
integer ldb,
real tola,
real tolb,
integer k,
integer l,
real, dimension( ldu, * ) u,
integer ldu,
real, dimension( ldv, * ) v,
integer ldv,
real, dimension( ldq, * ) q,
integer ldq,
integer, dimension( * ) iwork,
real, dimension( * ) tau,
real, dimension( * ) work,
integer info )

SGGSVP

Download SGGSVP + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!>
!> This routine is deprecated and has been replaced by routine SGGSVP3.
!>
!> SGGSVP computes orthogonal matrices U, V and Q such that
!>
!>                    N-K-L  K    L
!>  U**T*A*Q =     K ( 0    A12  A13 )  if M-K-L >= 0;
!>                 L ( 0     0   A23 )
!>             M-K-L ( 0     0    0  )
!>
!>                  N-K-L  K    L
!>         =     K ( 0    A12  A13 )  if M-K-L < 0;
!>             M-K ( 0     0   A23 )
!>
!>                  N-K-L  K    L
!>  V**T*B*Q =   L ( 0     0   B13 )
!>             P-L ( 0     0    0  )
!>
!> where the K-by-K matrix A12 and L-by-L matrix B13 are nonsingular
!> upper triangular; A23 is L-by-L upper triangular if M-K-L >= 0,
!> otherwise A23 is (M-K)-by-L upper trapezoidal.  K+L = the effective
!> numerical rank of the (M+P)-by-N matrix (A**T,B**T)**T.
!>
!> This decomposition is the preprocessing step for computing the
!> Generalized Singular Value Decomposition (GSVD), see subroutine
!> SGGSVD.
!> 
Parameters
[in]JOBU
!>          JOBU is CHARACTER*1
!>          = 'U':  Orthogonal matrix U is computed;
!>          = 'N':  U is not computed.
!> 
[in]JOBV
!>          JOBV is CHARACTER*1
!>          = 'V':  Orthogonal matrix V is computed;
!>          = 'N':  V is not computed.
!> 
[in]JOBQ
!>          JOBQ is CHARACTER*1
!>          = 'Q':  Orthogonal matrix Q is computed;
!>          = 'N':  Q is not computed.
!> 
[in]M
!>          M is INTEGER
!>          The number of rows of the matrix A.  M >= 0.
!> 
[in]P
!>          P is INTEGER
!>          The number of rows of the matrix B.  P >= 0.
!> 
[in]N
!>          N is INTEGER
!>          The number of columns of the matrices A and B.  N >= 0.
!> 
[in,out]A
!>          A is REAL array, dimension (LDA,N)
!>          On entry, the M-by-N matrix A.
!>          On exit, A contains the triangular (or trapezoidal) matrix
!>          described in the Purpose section.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A. LDA >= max(1,M).
!> 
[in,out]B
!>          B is REAL array, dimension (LDB,N)
!>          On entry, the P-by-N matrix B.
!>          On exit, B contains the triangular matrix described in
!>          the Purpose section.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B. LDB >= max(1,P).
!> 
[in]TOLA
!>          TOLA is REAL
!> 
[in]TOLB
!>          TOLB is REAL
!>
!>          TOLA and TOLB are the thresholds to determine the effective
!>          numerical rank of matrix B and a subblock of A. Generally,
!>          they are set to
!>             TOLA = MAX(M,N)*norm(A)*MACHEPS,
!>             TOLB = MAX(P,N)*norm(B)*MACHEPS.
!>          The size of TOLA and TOLB may affect the size of backward
!>          errors of the decomposition.
!> 
[out]K
!>          K is INTEGER
!> 
[out]L
!>          L is INTEGER
!>
!>          On exit, K and L specify the dimension of the subblocks
!>          described in Purpose section.
!>          K + L = effective numerical rank of (A**T,B**T)**T.
!> 
[out]U
!>          U is REAL array, dimension (LDU,M)
!>          If JOBU = 'U', U contains the orthogonal matrix U.
!>          If JOBU = 'N', U is not referenced.
!> 
[in]LDU
!>          LDU is INTEGER
!>          The leading dimension of the array U. LDU >= max(1,M) if
!>          JOBU = 'U'; LDU >= 1 otherwise.
!> 
[out]V
!>          V is REAL array, dimension (LDV,P)
!>          If JOBV = 'V', V contains the orthogonal matrix V.
!>          If JOBV = 'N', V is not referenced.
!> 
[in]LDV
!>          LDV is INTEGER
!>          The leading dimension of the array V. LDV >= max(1,P) if
!>          JOBV = 'V'; LDV >= 1 otherwise.
!> 
[out]Q
!>          Q is REAL array, dimension (LDQ,N)
!>          If JOBQ = 'Q', Q contains the orthogonal matrix Q.
!>          If JOBQ = 'N', Q is not referenced.
!> 
[in]LDQ
!>          LDQ is INTEGER
!>          The leading dimension of the array Q. LDQ >= max(1,N) if
!>          JOBQ = 'Q'; LDQ >= 1 otherwise.
!> 
[out]IWORK
!>          IWORK is INTEGER array, dimension (N)
!> 
[out]TAU
!>          TAU is REAL array, dimension (N)
!> 
[out]WORK
!>          WORK is REAL array, dimension (max(3*N,M,P))
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit
!>          < 0:  if INFO = -i, the i-th argument had an illegal value.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
The subroutine uses LAPACK subroutine SGEQPF for the QR factorization with column pivoting to detect the effective numerical rank of the a matrix. It may be replaced by a better rank determination strategy.

Definition at line 251 of file sggsvp.f.

254*
255* -- LAPACK computational routine --
256* -- LAPACK is a software package provided by Univ. of Tennessee, --
257* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
258*
259* .. Scalar Arguments ..
260 CHARACTER JOBQ, JOBU, JOBV
261 INTEGER INFO, K, L, LDA, LDB, LDQ, LDU, LDV, M, N, P
262 REAL TOLA, TOLB
263* ..
264* .. Array Arguments ..
265 INTEGER IWORK( * )
266 REAL A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
267 $ TAU( * ), U( LDU, * ), V( LDV, * ), WORK( * )
268* ..
269*
270* =====================================================================
271*
272* .. Parameters ..
273 REAL ZERO, ONE
274 parameter( zero = 0.0e+0, one = 1.0e+0 )
275* ..
276* .. Local Scalars ..
277 LOGICAL FORWRD, WANTQ, WANTU, WANTV
278 INTEGER I, J
279* ..
280* .. External Functions ..
281 LOGICAL LSAME
282 EXTERNAL lsame
283* ..
284* .. External Subroutines ..
285 EXTERNAL sgeqpf, sgeqr2, sgerq2, slacpy, slapmt, slaset,
287* ..
288* .. Intrinsic Functions ..
289 INTRINSIC abs, max, min
290* ..
291* .. Executable Statements ..
292*
293* Test the input parameters
294*
295 wantu = lsame( jobu, 'U' )
296 wantv = lsame( jobv, 'V' )
297 wantq = lsame( jobq, 'Q' )
298 forwrd = .true.
299*
300 info = 0
301 IF( .NOT.( wantu .OR. lsame( jobu, 'N' ) ) ) THEN
302 info = -1
303 ELSE IF( .NOT.( wantv .OR. lsame( jobv, 'N' ) ) ) THEN
304 info = -2
305 ELSE IF( .NOT.( wantq .OR. lsame( jobq, 'N' ) ) ) THEN
306 info = -3
307 ELSE IF( m.LT.0 ) THEN
308 info = -4
309 ELSE IF( p.LT.0 ) THEN
310 info = -5
311 ELSE IF( n.LT.0 ) THEN
312 info = -6
313 ELSE IF( lda.LT.max( 1, m ) ) THEN
314 info = -8
315 ELSE IF( ldb.LT.max( 1, p ) ) THEN
316 info = -10
317 ELSE IF( ldu.LT.1 .OR. ( wantu .AND. ldu.LT.m ) ) THEN
318 info = -16
319 ELSE IF( ldv.LT.1 .OR. ( wantv .AND. ldv.LT.p ) ) THEN
320 info = -18
321 ELSE IF( ldq.LT.1 .OR. ( wantq .AND. ldq.LT.n ) ) THEN
322 info = -20
323 END IF
324 IF( info.NE.0 ) THEN
325 CALL xerbla( 'SGGSVP', -info )
326 RETURN
327 END IF
328*
329* QR with column pivoting of B: B*P = V*( S11 S12 )
330* ( 0 0 )
331*
332 DO 10 i = 1, n
333 iwork( i ) = 0
334 10 CONTINUE
335 CALL sgeqpf( p, n, b, ldb, iwork, tau, work, info )
336*
337* Update A := A*P
338*
339 CALL slapmt( forwrd, m, n, a, lda, iwork )
340*
341* Determine the effective rank of matrix B.
342*
343 l = 0
344 DO 20 i = 1, min( p, n )
345 IF( abs( b( i, i ) ).GT.tolb )
346 $ l = l + 1
347 20 CONTINUE
348*
349 IF( wantv ) THEN
350*
351* Copy the details of V, and form V.
352*
353 CALL slaset( 'Full', p, p, zero, zero, v, ldv )
354 IF( p.GT.1 )
355 $ CALL slacpy( 'Lower', p-1, n, b( 2, 1 ), ldb, v( 2, 1 ),
356 $ ldv )
357 CALL sorg2r( p, p, min( p, n ), v, ldv, tau, work, info )
358 END IF
359*
360* Clean up B
361*
362 DO 40 j = 1, l - 1
363 DO 30 i = j + 1, l
364 b( i, j ) = zero
365 30 CONTINUE
366 40 CONTINUE
367 IF( p.GT.l )
368 $ CALL slaset( 'Full', p-l, n, zero, zero, b( l+1, 1 ), ldb )
369*
370 IF( wantq ) THEN
371*
372* Set Q = I and Update Q := Q*P
373*
374 CALL slaset( 'Full', n, n, zero, one, q, ldq )
375 CALL slapmt( forwrd, n, n, q, ldq, iwork )
376 END IF
377*
378 IF( p.GE.l .AND. n.NE.l ) THEN
379*
380* RQ factorization of (S11 S12): ( S11 S12 ) = ( 0 S12 )*Z
381*
382 CALL sgerq2( l, n, b, ldb, tau, work, info )
383*
384* Update A := A*Z**T
385*
386 CALL sormr2( 'Right', 'Transpose', m, n, l, b, ldb, tau, a,
387 $ lda, work, info )
388*
389 IF( wantq ) THEN
390*
391* Update Q := Q*Z**T
392*
393 CALL sormr2( 'Right', 'Transpose', n, n, l, b, ldb, tau,
394 $ q, ldq, work, info )
395 END IF
396*
397* Clean up B
398*
399 CALL slaset( 'Full', l, n-l, zero, zero, b, ldb )
400 DO 60 j = n - l + 1, n
401 DO 50 i = j - n + l + 1, l
402 b( i, j ) = zero
403 50 CONTINUE
404 60 CONTINUE
405*
406 END IF
407*
408* Let N-L L
409* A = ( A11 A12 ) M,
410*
411* then the following does the complete QR decomposition of A11:
412*
413* A11 = U*( 0 T12 )*P1**T
414* ( 0 0 )
415*
416 DO 70 i = 1, n - l
417 iwork( i ) = 0
418 70 CONTINUE
419 CALL sgeqpf( m, n-l, a, lda, iwork, tau, work, info )
420*
421* Determine the effective rank of A11
422*
423 k = 0
424 DO 80 i = 1, min( m, n-l )
425 IF( abs( a( i, i ) ).GT.tola )
426 $ k = k + 1
427 80 CONTINUE
428*
429* Update A12 := U**T*A12, where A12 = A( 1:M, N-L+1:N )
430*
431 CALL sorm2r( 'Left', 'Transpose', m, l, min( m, n-l ), a, lda,
432 $ tau, a( 1, n-l+1 ), lda, work, info )
433*
434 IF( wantu ) THEN
435*
436* Copy the details of U, and form U
437*
438 CALL slaset( 'Full', m, m, zero, zero, u, ldu )
439 IF( m.GT.1 )
440 $ CALL slacpy( 'Lower', m-1, n-l, a( 2, 1 ), lda,
441 $ u( 2, 1 ), ldu )
442 CALL sorg2r( m, m, min( m, n-l ), u, ldu, tau, work, info )
443 END IF
444*
445 IF( wantq ) THEN
446*
447* Update Q( 1:N, 1:N-L ) = Q( 1:N, 1:N-L )*P1
448*
449 CALL slapmt( forwrd, n, n-l, q, ldq, iwork )
450 END IF
451*
452* Clean up A: set the strictly lower triangular part of
453* A(1:K, 1:K) = 0, and A( K+1:M, 1:N-L ) = 0.
454*
455 DO 100 j = 1, k - 1
456 DO 90 i = j + 1, k
457 a( i, j ) = zero
458 90 CONTINUE
459 100 CONTINUE
460 IF( m.GT.k )
461 $ CALL slaset( 'Full', m-k, n-l, zero, zero, a( k+1, 1 ),
462 $ lda )
463*
464 IF( n-l.GT.k ) THEN
465*
466* RQ factorization of ( T11 T12 ) = ( 0 T12 )*Z1
467*
468 CALL sgerq2( k, n-l, a, lda, tau, work, info )
469*
470 IF( wantq ) THEN
471*
472* Update Q( 1:N,1:N-L ) = Q( 1:N,1:N-L )*Z1**T
473*
474 CALL sormr2( 'Right', 'Transpose', n, n-l, k, a, lda,
475 $ tau, q, ldq, work, info )
476 END IF
477*
478* Clean up A
479*
480 CALL slaset( 'Full', k, n-l-k, zero, zero, a, lda )
481 DO 120 j = n - l - k + 1, n - l
482 DO 110 i = j - n + l + k + 1, k
483 a( i, j ) = zero
484 110 CONTINUE
485 120 CONTINUE
486*
487 END IF
488*
489 IF( m.GT.k ) THEN
490*
491* QR factorization of A( K+1:M,N-L+1:N )
492*
493 CALL sgeqr2( m-k, l, a( k+1, n-l+1 ), lda, tau, work, info )
494*
495 IF( wantu ) THEN
496*
497* Update U(:,K+1:M) := U(:,K+1:M)*U1
498*
499 CALL sorm2r( 'Right', 'No transpose', m, m-k,
500 $ min( m-k, l ), a( k+1, n-l+1 ), lda, tau,
501 $ u( 1, k+1 ), ldu, work, info )
502 END IF
503*
504* Clean up
505*
506 DO 140 j = n - l + 1, n
507 DO 130 i = j - n + k + l + 1, m
508 a( i, j ) = zero
509 130 CONTINUE
510 140 CONTINUE
511*
512 END IF
513*
514 RETURN
515*
516* End of SGGSVP
517*
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:128
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:121
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:101
subroutine slapmt(forwrd, m, n, x, ldx, k)
SLAPMT performs a forward or backward permutation of the columns of a matrix.
Definition slapmt.f:102
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:108
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
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:112
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:157
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:157
subroutine sgeqpf(m, n, a, lda, jpvt, tau, work, info)
SGEQPF
Definition sgeqpf.f:140
Here is the call graph for this function:
Here is the caller graph for this function: