LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine cggsvp ( character  JOBU,
character  JOBV,
character  JOBQ,
integer  M,
integer  P,
integer  N,
complex, dimension( lda, * )  A,
integer  LDA,
complex, dimension( ldb, * )  B,
integer  LDB,
real  TOLA,
real  TOLB,
integer  K,
integer  L,
complex, dimension( ldu, * )  U,
integer  LDU,
complex, dimension( ldv, * )  V,
integer  LDV,
complex, dimension( ldq, * )  Q,
integer  LDQ,
integer, dimension( * )  IWORK,
real, dimension( * )  RWORK,
complex, dimension( * )  TAU,
complex, dimension( * )  WORK,
integer  INFO 
)

CGGSVP

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

Purpose:
 This routine is deprecated and has been replaced by routine CGGSVP3.

 CGGSVP computes unitary matrices U, V and Q such that

                    N-K-L  K    L
  U**H*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**H*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**H,B**H)**H. 

 This decomposition is the preprocessing step for computing the
 Generalized Singular Value Decomposition (GSVD), see subroutine
 CGGSVD.
Parameters
[in]JOBU
          JOBU is CHARACTER*1
          = 'U':  Unitary matrix U is computed;
          = 'N':  U is not computed.
[in]JOBV
          JOBV is CHARACTER*1
          = 'V':  Unitary matrix V is computed;
          = 'N':  V is not computed.
[in]JOBQ
          JOBQ is CHARACTER*1
          = 'Q':  Unitary 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 COMPLEX 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 COMPLEX 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**H,B**H)**H.
[out]U
          U is COMPLEX array, dimension (LDU,M)
          If JOBU = 'U', U contains the unitary 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 COMPLEX array, dimension (LDV,P)
          If JOBV = 'V', V contains the unitary 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 COMPLEX array, dimension (LDQ,N)
          If JOBQ = 'Q', Q contains the unitary 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]RWORK
          RWORK is REAL array, dimension (2*N)
[out]TAU
          TAU is COMPLEX array, dimension (N)
[out]WORK
          WORK is COMPLEX 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.
Date
November 2011
Further Details:
The subroutine uses LAPACK subroutine CGEQPF 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 264 of file cggsvp.f.

264 *
265 * -- LAPACK computational routine (version 3.4.0) --
266 * -- LAPACK is a software package provided by Univ. of Tennessee, --
267 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
268 * November 2011
269 *
270 * .. Scalar Arguments ..
271  CHARACTER jobq, jobu, jobv
272  INTEGER info, k, l, lda, ldb, ldq, ldu, ldv, m, n, p
273  REAL tola, tolb
274 * ..
275 * .. Array Arguments ..
276  INTEGER iwork( * )
277  REAL rwork( * )
278  COMPLEX a( lda, * ), b( ldb, * ), q( ldq, * ),
279  $ tau( * ), u( ldu, * ), v( ldv, * ), work( * )
280 * ..
281 *
282 * =====================================================================
283 *
284 * .. Parameters ..
285  COMPLEX czero, cone
286  parameter ( czero = ( 0.0e+0, 0.0e+0 ),
287  $ cone = ( 1.0e+0, 0.0e+0 ) )
288 * ..
289 * .. Local Scalars ..
290  LOGICAL forwrd, wantq, wantu, wantv
291  INTEGER i, j
292  COMPLEX t
293 * ..
294 * .. External Functions ..
295  LOGICAL lsame
296  EXTERNAL lsame
297 * ..
298 * .. External Subroutines ..
299  EXTERNAL cgeqpf, cgeqr2, cgerq2, clacpy, clapmt, claset,
301 * ..
302 * .. Intrinsic Functions ..
303  INTRINSIC abs, aimag, max, min, real
304 * ..
305 * .. Statement Functions ..
306  REAL cabs1
307 * ..
308 * .. Statement Function definitions ..
309  cabs1( t ) = abs( REAL( T ) ) + abs( aimag( t ) )
310 * ..
311 * .. Executable Statements ..
312 *
313 * Test the input parameters
314 *
315  wantu = lsame( jobu, 'U' )
316  wantv = lsame( jobv, 'V' )
317  wantq = lsame( jobq, 'Q' )
318  forwrd = .true.
319 *
320  info = 0
321  IF( .NOT.( wantu .OR. lsame( jobu, 'N' ) ) ) THEN
322  info = -1
323  ELSE IF( .NOT.( wantv .OR. lsame( jobv, 'N' ) ) ) THEN
324  info = -2
325  ELSE IF( .NOT.( wantq .OR. lsame( jobq, 'N' ) ) ) THEN
326  info = -3
327  ELSE IF( m.LT.0 ) THEN
328  info = -4
329  ELSE IF( p.LT.0 ) THEN
330  info = -5
331  ELSE IF( n.LT.0 ) THEN
332  info = -6
333  ELSE IF( lda.LT.max( 1, m ) ) THEN
334  info = -8
335  ELSE IF( ldb.LT.max( 1, p ) ) THEN
336  info = -10
337  ELSE IF( ldu.LT.1 .OR. ( wantu .AND. ldu.LT.m ) ) THEN
338  info = -16
339  ELSE IF( ldv.LT.1 .OR. ( wantv .AND. ldv.LT.p ) ) THEN
340  info = -18
341  ELSE IF( ldq.LT.1 .OR. ( wantq .AND. ldq.LT.n ) ) THEN
342  info = -20
343  END IF
344  IF( info.NE.0 ) THEN
345  CALL xerbla( 'CGGSVP', -info )
346  RETURN
347  END IF
348 *
349 * QR with column pivoting of B: B*P = V*( S11 S12 )
350 * ( 0 0 )
351 *
352  DO 10 i = 1, n
353  iwork( i ) = 0
354  10 CONTINUE
355  CALL cgeqpf( p, n, b, ldb, iwork, tau, work, rwork, info )
356 *
357 * Update A := A*P
358 *
359  CALL clapmt( forwrd, m, n, a, lda, iwork )
360 *
361 * Determine the effective rank of matrix B.
362 *
363  l = 0
364  DO 20 i = 1, min( p, n )
365  IF( cabs1( b( i, i ) ).GT.tolb )
366  $ l = l + 1
367  20 CONTINUE
368 *
369  IF( wantv ) THEN
370 *
371 * Copy the details of V, and form V.
372 *
373  CALL claset( 'Full', p, p, czero, czero, v, ldv )
374  IF( p.GT.1 )
375  $ CALL clacpy( 'Lower', p-1, n, b( 2, 1 ), ldb, v( 2, 1 ),
376  $ ldv )
377  CALL cung2r( p, p, min( p, n ), v, ldv, tau, work, info )
378  END IF
379 *
380 * Clean up B
381 *
382  DO 40 j = 1, l - 1
383  DO 30 i = j + 1, l
384  b( i, j ) = czero
385  30 CONTINUE
386  40 CONTINUE
387  IF( p.GT.l )
388  $ CALL claset( 'Full', p-l, n, czero, czero, b( l+1, 1 ), ldb )
389 *
390  IF( wantq ) THEN
391 *
392 * Set Q = I and Update Q := Q*P
393 *
394  CALL claset( 'Full', n, n, czero, cone, q, ldq )
395  CALL clapmt( forwrd, n, n, q, ldq, iwork )
396  END IF
397 *
398  IF( p.GE.l .AND. n.NE.l ) THEN
399 *
400 * RQ factorization of ( S11 S12 ) = ( 0 S12 )*Z
401 *
402  CALL cgerq2( l, n, b, ldb, tau, work, info )
403 *
404 * Update A := A*Z**H
405 *
406  CALL cunmr2( 'Right', 'Conjugate transpose', m, n, l, b, ldb,
407  $ tau, a, lda, work, info )
408  IF( wantq ) THEN
409 *
410 * Update Q := Q*Z**H
411 *
412  CALL cunmr2( 'Right', 'Conjugate transpose', n, n, l, b,
413  $ ldb, tau, q, ldq, work, info )
414  END IF
415 *
416 * Clean up B
417 *
418  CALL claset( 'Full', l, n-l, czero, czero, b, ldb )
419  DO 60 j = n - l + 1, n
420  DO 50 i = j - n + l + 1, l
421  b( i, j ) = czero
422  50 CONTINUE
423  60 CONTINUE
424 *
425  END IF
426 *
427 * Let N-L L
428 * A = ( A11 A12 ) M,
429 *
430 * then the following does the complete QR decomposition of A11:
431 *
432 * A11 = U*( 0 T12 )*P1**H
433 * ( 0 0 )
434 *
435  DO 70 i = 1, n - l
436  iwork( i ) = 0
437  70 CONTINUE
438  CALL cgeqpf( m, n-l, a, lda, iwork, tau, work, rwork, info )
439 *
440 * Determine the effective rank of A11
441 *
442  k = 0
443  DO 80 i = 1, min( m, n-l )
444  IF( cabs1( a( i, i ) ).GT.tola )
445  $ k = k + 1
446  80 CONTINUE
447 *
448 * Update A12 := U**H*A12, where A12 = A( 1:M, N-L+1:N )
449 *
450  CALL cunm2r( 'Left', 'Conjugate transpose', m, l, min( m, n-l ),
451  $ a, lda, tau, a( 1, n-l+1 ), lda, work, info )
452 *
453  IF( wantu ) THEN
454 *
455 * Copy the details of U, and form U
456 *
457  CALL claset( 'Full', m, m, czero, czero, u, ldu )
458  IF( m.GT.1 )
459  $ CALL clacpy( 'Lower', m-1, n-l, a( 2, 1 ), lda, u( 2, 1 ),
460  $ ldu )
461  CALL cung2r( m, m, min( m, n-l ), u, ldu, tau, work, info )
462  END IF
463 *
464  IF( wantq ) THEN
465 *
466 * Update Q( 1:N, 1:N-L ) = Q( 1:N, 1:N-L )*P1
467 *
468  CALL clapmt( forwrd, n, n-l, q, ldq, iwork )
469  END IF
470 *
471 * Clean up A: set the strictly lower triangular part of
472 * A(1:K, 1:K) = 0, and A( K+1:M, 1:N-L ) = 0.
473 *
474  DO 100 j = 1, k - 1
475  DO 90 i = j + 1, k
476  a( i, j ) = czero
477  90 CONTINUE
478  100 CONTINUE
479  IF( m.GT.k )
480  $ CALL claset( 'Full', m-k, n-l, czero, czero, a( k+1, 1 ), lda )
481 *
482  IF( n-l.GT.k ) THEN
483 *
484 * RQ factorization of ( T11 T12 ) = ( 0 T12 )*Z1
485 *
486  CALL cgerq2( k, n-l, a, lda, tau, work, info )
487 *
488  IF( wantq ) THEN
489 *
490 * Update Q( 1:N,1:N-L ) = Q( 1:N,1:N-L )*Z1**H
491 *
492  CALL cunmr2( 'Right', 'Conjugate transpose', n, n-l, k, a,
493  $ lda, tau, q, ldq, work, info )
494  END IF
495 *
496 * Clean up A
497 *
498  CALL claset( 'Full', k, n-l-k, czero, czero, a, lda )
499  DO 120 j = n - l - k + 1, n - l
500  DO 110 i = j - n + l + k + 1, k
501  a( i, j ) = czero
502  110 CONTINUE
503  120 CONTINUE
504 *
505  END IF
506 *
507  IF( m.GT.k ) THEN
508 *
509 * QR factorization of A( K+1:M,N-L+1:N )
510 *
511  CALL cgeqr2( m-k, l, a( k+1, n-l+1 ), lda, tau, work, info )
512 *
513  IF( wantu ) THEN
514 *
515 * Update U(:,K+1:M) := U(:,K+1:M)*U1
516 *
517  CALL cunm2r( 'Right', 'No transpose', m, m-k, min( m-k, l ),
518  $ a( k+1, n-l+1 ), lda, tau, u( 1, k+1 ), ldu,
519  $ work, info )
520  END IF
521 *
522 * Clean up
523 *
524  DO 140 j = n - l + 1, n
525  DO 130 i = j - n + k + l + 1, m
526  a( i, j ) = czero
527  130 CONTINUE
528  140 CONTINUE
529 *
530  END IF
531 *
532  RETURN
533 *
534 * End of CGGSVP
535 *
subroutine cgeqr2(M, N, A, LDA, TAU, WORK, INFO)
CGEQR2 computes the QR factorization of a general rectangular matrix using an unblocked algorithm...
Definition: cgeqr2.f:123
subroutine cunm2r(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, INFO)
CUNM2R multiplies a general matrix by the unitary matrix from a QR factorization determined by cgeqrf...
Definition: cunm2r.f:161
subroutine cung2r(M, N, K, A, LDA, TAU, WORK, INFO)
CUNG2R
Definition: cung2r.f:116
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine cgeqpf(M, N, A, LDA, JPVT, TAU, WORK, RWORK, INFO)
CGEQPF
Definition: cgeqpf.f:150
subroutine cunmr2(SIDE, TRANS, M, N, K, A, LDA, TAU, C, LDC, WORK, INFO)
CUNMR2 multiplies a general matrix by the unitary matrix from a RQ factorization determined by cgerqf...
Definition: cunmr2.f:161
subroutine claset(UPLO, M, N, ALPHA, BETA, A, LDA)
CLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: claset.f:108
subroutine clacpy(UPLO, M, N, A, LDA, B, LDB)
CLACPY copies all or part of one two-dimensional array to another.
Definition: clacpy.f:105
subroutine cgerq2(M, N, A, LDA, TAU, WORK, INFO)
CGERQ2 computes the RQ factorization of a general rectangular matrix using an unblocked algorithm...
Definition: cgerq2.f:125
subroutine clapmt(FORWRD, M, N, X, LDX, K)
CLAPMT performs a forward or backward permutation of the columns of a matrix.
Definition: clapmt.f:106
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: