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

◆ cggsvp3()

subroutine cggsvp3 ( 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 lwork,
integer info )

CGGSVP3

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

Purpose:
!>
!> CGGSVP3 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
!> CGGSVD3.
!> 
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(1,LWORK))
!>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
!> 
[in]LWORK
!>          LWORK is INTEGER
!>          The dimension of the array WORK. LWORK >= 1.
!>
!>          If LWORK = -1, then a workspace query is assumed; the routine
!>          only calculates the optimal size of the WORK array, returns
!>          this value as the first entry of the WORK array, and no error
!>          message related to LWORK is issued by XERBLA.
!> 
[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 CGEQP3 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.
!>
!>  CGGSVP3 replaces the deprecated subroutine CGGSVP.
!>
!> 

Definition at line 273 of file cggsvp3.f.

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