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

◆ cggsvp()

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.
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 259 of file cggsvp.f.

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