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

◆ zggsvd3()

subroutine zggsvd3 ( character jobu,
character jobv,
character jobq,
integer m,
integer n,
integer p,
integer k,
integer l,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( ldb, * ) b,
integer ldb,
double precision, dimension( * ) alpha,
double precision, dimension( * ) beta,
complex*16, dimension( ldu, * ) u,
integer ldu,
complex*16, dimension( ldv, * ) v,
integer ldv,
complex*16, dimension( ldq, * ) q,
integer ldq,
complex*16, dimension( * ) work,
integer lwork,
double precision, dimension( * ) rwork,
integer, dimension( * ) iwork,
integer info )

ZGGSVD3 computes the singular value decomposition (SVD) for OTHER matrices

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

Purpose:
!>
!> ZGGSVD3 computes the generalized singular value decomposition (GSVD)
!> of an M-by-N complex matrix A and P-by-N complex matrix B:
!>
!>       U**H*A*Q = D1*( 0 R ),    V**H*B*Q = D2*( 0 R )
!>
!> where U, V and Q are unitary matrices.
!> Let K+L = the effective numerical rank of the
!> matrix (A**H,B**H)**H, then R is a (K+L)-by-(K+L) nonsingular upper
!> triangular matrix, D1 and D2 are M-by-(K+L) and P-by-(K+L) 
!> matrices and of the following structures, respectively:
!>
!> If M-K-L >= 0,
!>
!>                     K  L
!>        D1 =     K ( I  0 )
!>                 L ( 0  C )
!>             M-K-L ( 0  0 )
!>
!>                   K  L
!>        D2 =   L ( 0  S )
!>             P-L ( 0  0 )
!>
!>                 N-K-L  K    L
!>   ( 0 R ) = K (  0   R11  R12 )
!>             L (  0    0   R22 )
!> where
!>
!>   C = diag( ALPHA(K+1), ... , ALPHA(K+L) ),
!>   S = diag( BETA(K+1),  ... , BETA(K+L) ),
!>   C**2 + S**2 = I.
!>
!>   R is stored in A(1:K+L,N-K-L+1:N) on exit.
!>
!> If M-K-L < 0,
!>
!>                   K M-K K+L-M
!>        D1 =   K ( I  0    0   )
!>             M-K ( 0  C    0   )
!>
!>                     K M-K K+L-M
!>        D2 =   M-K ( 0  S    0  )
!>             K+L-M ( 0  0    I  )
!>               P-L ( 0  0    0  )
!>
!>                    N-K-L  K   M-K  K+L-M
!>   ( 0 R ) =     K ( 0    R11  R12  R13  )
!>               M-K ( 0     0   R22  R23  )
!>             K+L-M ( 0     0    0   R33  )
!>
!> where
!>
!>   C = diag( ALPHA(K+1), ... , ALPHA(M) ),
!>   S = diag( BETA(K+1),  ... , BETA(M) ),
!>   C**2 + S**2 = I.
!>
!>   (R11 R12 R13 ) is stored in A(1:M, N-K-L+1:N), and R33 is stored
!>   ( 0  R22 R23 )
!>   in B(M-K+1:L,N+M-K-L+1:N) on exit.
!>
!> The routine computes C, S, R, and optionally the unitary
!> transformation matrices U, V and Q.
!>
!> In particular, if B is an N-by-N nonsingular matrix, then the GSVD of
!> A and B implicitly gives the SVD of A*inv(B):
!>                      A*inv(B) = U*(D1*inv(D2))*V**H.
!> If ( A**H,B**H)**H has orthonormal columns, then the GSVD of A and B is also
!> equal to the CS decomposition of A and B. Furthermore, the GSVD can
!> be used to derive the solution of the eigenvalue problem:
!>                      A**H*A x = lambda* B**H*B x.
!> In some literature, the GSVD of A and B is presented in the form
!>                  U**H*A*X = ( 0 D1 ),   V**H*B*X = ( 0 D2 )
!> where U and V are orthogonal and X is nonsingular, and D1 and D2 are
!> ``diagonal''.  The former GSVD form can be converted to the latter
!> form by taking the nonsingular matrix X as
!>
!>                       X = Q*(  I   0    )
!>                             (  0 inv(R) )
!> 
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]N
!>          N is INTEGER
!>          The number of columns of the matrices A and B.  N >= 0.
!> 
[in]P
!>          P is INTEGER
!>          The number of rows of the matrix B.  P >= 0.
!> 
[out]K
!>          K is INTEGER
!> 
[out]L
!>          L is INTEGER
!>
!>          On exit, K and L specify the dimension of the subblocks
!>          described in Purpose.
!>          K + L = effective numerical rank of (A**H,B**H)**H.
!> 
[in,out]A
!>          A is COMPLEX*16 array, dimension (LDA,N)
!>          On entry, the M-by-N matrix A.
!>          On exit, A contains the triangular matrix R, or part of R.
!>          See Purpose for details.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A. LDA >= max(1,M).
!> 
[in,out]B
!>          B is COMPLEX*16 array, dimension (LDB,N)
!>          On entry, the P-by-N matrix B.
!>          On exit, B contains part of the triangular matrix R if
!>          M-K-L < 0.  See Purpose for details.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B. LDB >= max(1,P).
!> 
[out]ALPHA
!>          ALPHA is DOUBLE PRECISION array, dimension (N)
!> 
[out]BETA
!>          BETA is DOUBLE PRECISION array, dimension (N)
!>
!>          On exit, ALPHA and BETA contain the generalized singular
!>          value pairs of A and B;
!>            ALPHA(1:K) = 1,
!>            BETA(1:K)  = 0,
!>          and if M-K-L >= 0,
!>            ALPHA(K+1:K+L) = C,
!>            BETA(K+1:K+L)  = S,
!>          or if M-K-L < 0,
!>            ALPHA(K+1:M)=C, ALPHA(M+1:K+L)=0
!>            BETA(K+1:M) =S, BETA(M+1:K+L) =1
!>          and
!>            ALPHA(K+L+1:N) = 0
!>            BETA(K+L+1:N)  = 0
!> 
[out]U
!>          U is COMPLEX*16 array, dimension (LDU,M)
!>          If JOBU = 'U', U contains the M-by-M 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*16 array, dimension (LDV,P)
!>          If JOBV = 'V', V contains the P-by-P 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*16 array, dimension (LDQ,N)
!>          If JOBQ = 'Q', Q contains the N-by-N 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]WORK
!>          WORK is COMPLEX*16 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]RWORK
!>          RWORK is DOUBLE PRECISION array, dimension (2*N)
!> 
[out]IWORK
!>          IWORK is INTEGER array, dimension (N)
!>          On exit, IWORK stores the sorting information. More
!>          precisely, the following loop will sort ALPHA
!>             for I = K+1, min(M,K+L)
!>                 swap ALPHA(I) and ALPHA(IWORK(I))
!>             endfor
!>          such that ALPHA(1) >= ALPHA(2) >= ... >= ALPHA(N).
!> 
[out]INFO
!>          INFO is INTEGER
!>          = 0:  successful exit.
!>          < 0:  if INFO = -i, the i-th argument had an illegal value.
!>          > 0:  if INFO = 1, the Jacobi-type procedure failed to
!>                converge.  For further details, see subroutine ZTGSJA.
!> 
Internal Parameters:
!>  TOLA    DOUBLE PRECISION
!>  TOLB    DOUBLE PRECISION
!>          TOLA and TOLB are the thresholds to determine the effective
!>          rank of (A**H,B**H)**H. 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.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Ming Gu and Huan Ren, Computer Science Division, University of California at Berkeley, USA
Further Details:
ZGGSVD3 replaces the deprecated subroutine ZGGSVD.

Definition at line 348 of file zggsvd3.f.

351*
352* -- LAPACK driver routine --
353* -- LAPACK is a software package provided by Univ. of Tennessee, --
354* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
355*
356* .. Scalar Arguments ..
357 CHARACTER JOBQ, JOBU, JOBV
358 INTEGER INFO, K, L, LDA, LDB, LDQ, LDU, LDV, M, N, P,
359 $ LWORK
360* ..
361* .. Array Arguments ..
362 INTEGER IWORK( * )
363 DOUBLE PRECISION ALPHA( * ), BETA( * ), RWORK( * )
364 COMPLEX*16 A( LDA, * ), B( LDB, * ), Q( LDQ, * ),
365 $ U( LDU, * ), V( LDV, * ), WORK( * )
366* ..
367*
368* =====================================================================
369*
370* .. Local Scalars ..
371 LOGICAL WANTQ, WANTU, WANTV, LQUERY
372 INTEGER I, IBND, ISUB, J, NCYCLE, LWKOPT
373 DOUBLE PRECISION ANORM, BNORM, SMAX, TEMP, TOLA, TOLB, ULP, UNFL
374* ..
375* .. External Functions ..
376 LOGICAL LSAME
377 DOUBLE PRECISION DLAMCH, ZLANGE
378 EXTERNAL lsame, dlamch, zlange
379* ..
380* .. External Subroutines ..
381 EXTERNAL dcopy, xerbla, zggsvp3, ztgsja
382* ..
383* .. Intrinsic Functions ..
384 INTRINSIC max, min
385* ..
386* .. Executable Statements ..
387*
388* Decode and test the input parameters
389*
390 wantu = lsame( jobu, 'U' )
391 wantv = lsame( jobv, 'V' )
392 wantq = lsame( jobq, 'Q' )
393 lquery = ( lwork.EQ.-1 )
394 lwkopt = 1
395*
396* Test the input arguments
397*
398 info = 0
399 IF( .NOT.( wantu .OR. lsame( jobu, 'N' ) ) ) THEN
400 info = -1
401 ELSE IF( .NOT.( wantv .OR. lsame( jobv, 'N' ) ) ) THEN
402 info = -2
403 ELSE IF( .NOT.( wantq .OR. lsame( jobq, 'N' ) ) ) THEN
404 info = -3
405 ELSE IF( m.LT.0 ) THEN
406 info = -4
407 ELSE IF( n.LT.0 ) THEN
408 info = -5
409 ELSE IF( p.LT.0 ) THEN
410 info = -6
411 ELSE IF( lda.LT.max( 1, m ) ) THEN
412 info = -10
413 ELSE IF( ldb.LT.max( 1, p ) ) THEN
414 info = -12
415 ELSE IF( ldu.LT.1 .OR. ( wantu .AND. ldu.LT.m ) ) THEN
416 info = -16
417 ELSE IF( ldv.LT.1 .OR. ( wantv .AND. ldv.LT.p ) ) THEN
418 info = -18
419 ELSE IF( ldq.LT.1 .OR. ( wantq .AND. ldq.LT.n ) ) THEN
420 info = -20
421 ELSE IF( lwork.LT.1 .AND. .NOT.lquery ) THEN
422 info = -24
423 END IF
424*
425* Compute workspace
426*
427 IF( info.EQ.0 ) THEN
428 CALL zggsvp3( jobu, jobv, jobq, m, p, n, a, lda, b, ldb,
429 $ tola,
430 $ tolb, k, l, u, ldu, v, ldv, q, ldq, iwork, rwork,
431 $ work, work, -1, info )
432 lwkopt = n + int( work( 1 ) )
433 lwkopt = max( 2*n, lwkopt )
434 lwkopt = max( 1, lwkopt )
435 work( 1 ) = dcmplx( lwkopt )
436 END IF
437*
438 IF( info.NE.0 ) THEN
439 CALL xerbla( 'ZGGSVD3', -info )
440 RETURN
441 END IF
442 IF( lquery ) THEN
443 RETURN
444 ENDIF
445*
446* Compute the Frobenius norm of matrices A and B
447*
448 anorm = zlange( '1', m, n, a, lda, rwork )
449 bnorm = zlange( '1', p, n, b, ldb, rwork )
450*
451* Get machine precision and set up threshold for determining
452* the effective numerical rank of the matrices A and B.
453*
454 ulp = dlamch( 'Precision' )
455 unfl = dlamch( 'Safe Minimum' )
456 tola = max( m, n )*max( anorm, unfl )*ulp
457 tolb = max( p, n )*max( bnorm, unfl )*ulp
458*
459 CALL zggsvp3( jobu, jobv, jobq, m, p, n, a, lda, b, ldb, tola,
460 $ tolb, k, l, u, ldu, v, ldv, q, ldq, iwork, rwork,
461 $ work, work( n+1 ), lwork-n, info )
462*
463* Compute the GSVD of two upper "triangular" matrices
464*
465 CALL ztgsja( jobu, jobv, jobq, m, p, n, k, l, a, lda, b, ldb,
466 $ tola, tolb, alpha, beta, u, ldu, v, ldv, q, ldq,
467 $ work, ncycle, info )
468*
469* Sort the singular values and store the pivot indices in IWORK
470* Copy ALPHA to RWORK, then sort ALPHA in RWORK
471*
472 CALL dcopy( n, alpha, 1, rwork, 1 )
473 ibnd = min( l, m-k )
474 DO 20 i = 1, ibnd
475*
476* Scan for largest ALPHA(K+I)
477*
478 isub = i
479 smax = rwork( k+i )
480 DO 10 j = i + 1, ibnd
481 temp = rwork( k+j )
482 IF( temp.GT.smax ) THEN
483 isub = j
484 smax = temp
485 END IF
486 10 CONTINUE
487 IF( isub.NE.i ) THEN
488 rwork( k+isub ) = rwork( k+i )
489 rwork( k+i ) = smax
490 iwork( k+i ) = k + isub
491 ELSE
492 iwork( k+i ) = k + i
493 END IF
494 20 CONTINUE
495*
496 work( 1 ) = dcmplx( lwkopt )
497 RETURN
498*
499* End of ZGGSVD3
500*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
subroutine zggsvp3(jobu, jobv, jobq, m, p, n, a, lda, b, ldb, tola, tolb, k, l, u, ldu, v, ldv, q, ldq, iwork, rwork, tau, work, lwork, info)
ZGGSVP3
Definition zggsvp3.f:276
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlange(norm, m, n, a, lda, work)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlange.f:113
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine ztgsja(jobu, jobv, jobq, m, p, n, k, l, a, lda, b, ldb, tola, tolb, alpha, beta, u, ldu, v, ldv, q, ldq, work, ncycle, info)
ZTGSJA
Definition ztgsja.f:377
Here is the call graph for this function:
Here is the caller graph for this function: