NAG Library Routine Document

f08vaf (dggsvd)

1
Purpose

f08vaf (dggsvd) computes the generalized singular value decomposition (GSVD) of an m by n real matrix A and a p by n real matrix B. f08vaf (dggsvd) is marked as deprecated by LAPACK; the replacement routine is f08vcf (dggsvd3) which makes better use of level 3 BLAS.

2
Specification

Fortran Interface
Subroutine f08vaf ( jobu, jobv, jobq, m, n, p, k, l, a, lda, b, ldb, alpha, beta, u, ldu, v, ldv, q, ldq, work, iwork, info)
Integer, Intent (In):: m, n, p, lda, ldb, ldu, ldv, ldq
Integer, Intent (Out):: k, l, iwork(n), info
Real (Kind=nag_wp), Intent (Inout):: a(lda,*), b(ldb,*), u(ldu,*), v(ldv,*), q(ldq,*)
Real (Kind=nag_wp), Intent (Out):: alpha(n), beta(n), work(max(3*n,m,p)+n)
Character (1), Intent (In):: jobu, jobv, jobq
C Header Interface
#include <nagmk26.h>
void  f08vaf_ (const char *jobu, const char *jobv, const char *jobq, const Integer *m, const Integer *n, const Integer *p, Integer *k, Integer *l, double a[], const Integer *lda, double b[], const Integer *ldb, double alpha[], double beta[], double u[], const Integer *ldu, double v[], const Integer *ldv, double q[], const Integer *ldq, double work[], Integer iwork[], Integer *info, const Charlen length_jobu, const Charlen length_jobv, const Charlen length_jobq)
The routine may be called by its LAPACK name dggsvd.

3
Description

The generalized singular value decomposition is given by
UT A Q = D1 0 R ,   VT B Q = D2 0 R ,  
where U, V and Q are orthogonal matrices. Let k+l be the effective numerical rank of the matrix A B , 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 ‘diagonal’ matrices structured as follows:
if m-k-l0,
D1= klkI0l0Cm-k-l00()  
D2= kll0Sp-l00()  
0R = n-k-lklk0R11R12l00R22()  
where
C = diagαk+1,,αk+l ,  
S = diagβk+1,,βk+l ,  
and
C2 + S2 = I .  
R is stored as a submatrix of A with elements Rij stored as Ai,n-k-l+j on exit.
If m-k-l<0 ,
D1= km-kk+l-mkI00m-k0C0()  
D2= km-kk+l-mm-k0S0k+l-m00Ip-l000()  
0R = n-k-lkm-kk+l-mk0R11R12R13m-k00R22R23k+l-m000R33()  
where
C = diagαk+1,,αm ,  
S = diagβk+1,,βm ,  
and
C2 + S2 = I .  
R11 R12 R13 0 R22 R23  is stored as a submatrix of A with Rij stored as Ai,n-k-l+j, and R33  is stored as a submatrix of B with R33ij stored as Bm-k+i,n+m-k-l+j.
The routine computes C, S, R and, optionally, the orthogonal 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 AB-1:
A B-1 = U D1 D2-1 VT .  
If A B  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:
AT Ax=λ BT Bx .  
In some literature, the GSVD of A and B is presented in the form
UT A X = 0D1 ,   VT B X = 0D2 ,  
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 R-1 .  

4
References

Anderson E, Bai Z, Bischof C, Blackford S, Demmel J, Dongarra J J, Du Croz J J, Greenbaum A, Hammarling S, McKenney A and Sorensen D (1999) LAPACK Users' Guide (3rd Edition) SIAM, Philadelphia http://www.netlib.org/lapack/lug
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore

5
Arguments

1:     jobu – Character(1)Input
On entry: if jobu='U', the orthogonal matrix U is computed.
If jobu='N', U is not computed.
Constraint: jobu='U' or 'N'.
2:     jobv – Character(1)Input
On entry: if jobv='V', the orthogonal matrix V is computed.
If jobv='N', V is not computed.
Constraint: jobv='V' or 'N'.
3:     jobq – Character(1)Input
On entry: if jobq='Q', the orthogonal matrix Q is computed.
If jobq='N', Q is not computed.
Constraint: jobq='Q' or 'N'.
4:     m – IntegerInput
On entry: m, the number of rows of the matrix A.
Constraint: m0.
5:     n – IntegerInput
On entry: n, the number of columns of the matrices A and B.
Constraint: n0.
6:     p – IntegerInput
On entry: p, the number of rows of the matrix B.
Constraint: p0.
7:     k – IntegerOutput
8:     l – IntegerOutput
On exit: k and l specify the dimension of the subblocks k and l as described in Section 3; k+l is the effective numerical rank of AB.
9:     alda* – Real (Kind=nag_wp) arrayInput/Output
Note: the second dimension of the array a must be at least max1,n.
On entry: the m by n matrix A.
On exit: contains the triangular matrix R, or part of R. See Section 3 for details.
10:   lda – IntegerInput
On entry: the first dimension of the array a as declared in the (sub)program from which f08vaf (dggsvd) is called.
Constraint: ldamax1,m.
11:   bldb* – Real (Kind=nag_wp) arrayInput/Output
Note: the second dimension of the array b must be at least max1,n.
On entry: the p by n matrix B.
On exit: contains the triangular matrix R if m-k-l<0. See Section 3 for details.
12:   ldb – IntegerInput
On entry: the first dimension of the array b as declared in the (sub)program from which f08vaf (dggsvd) is called.
Constraint: ldbmax1,p.
13:   alphan – Real (Kind=nag_wp) arrayOutput
On exit: see the description of beta.
14:   betan – Real (Kind=nag_wp) arrayOutput
On exit: alpha and beta contain the generalized singular value pairs of A and B, αi  and βi ;
  • alpha1:k = 1 ,
  • beta1:k = 0 ,
and if m-k-l0 ,
  • alphak+1:k+l = C ,
  • betak+1:k+l = S ,
or if m-k-l<0 ,
  • alphak+1:m = C ,
  • alpham+1:k+l = 0 ,
  • betak+1:m = S ,
  • betam+1:k+l = 1 , and
  • alphak+l+1:n = 0 ,
  • betak+l+1:n = 0 .
The notation alphak:n above refers to consecutive elements alphai, for i=k,,n.
15:   uldu* – Real (Kind=nag_wp) arrayOutput
Note: the second dimension of the array u must be at least max1,m if jobu='U'.
On exit: if jobu='U', u contains the m by m orthogonal matrix U.
If jobu='N', u is not referenced.
16:   ldu – IntegerInput
On entry: the first dimension of the array u as declared in the (sub)program from which f08vaf (dggsvd) is called.
Constraints:
  • if jobu='U', ldu max1,m ;
  • otherwise ldu1.
17:   vldv* – Real (Kind=nag_wp) arrayOutput
Note: the second dimension of the array v must be at least max1,p if jobv='V'.
On exit: if jobv='V', v contains the p by p orthogonal matrix V.
If jobv='N', v is not referenced.
18:   ldv – IntegerInput
On entry: the first dimension of the array v as declared in the (sub)program from which f08vaf (dggsvd) is called.
Constraints:
  • if jobv='V', ldv max1,p ;
  • otherwise ldv1.
19:   qldq* – Real (Kind=nag_wp) arrayOutput
Note: the second dimension of the array q must be at least max1,n if jobq='Q'.
On exit: if jobq='Q', q contains the n by n orthogonal matrix Q.
If jobq='N', q is not referenced.
20:   ldq – IntegerInput
On entry: the first dimension of the array q as declared in the (sub)program from which f08vaf (dggsvd) is called.
Constraints:
  • if jobq='Q', ldq max1,n ;
  • otherwise ldq1.
21:   workmax3×n,m,p+n – Real (Kind=nag_wp) arrayWorkspace
22:   iworkn – Integer arrayOutput
On exit: 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 alpha1alpha2alphan.
23:   info – IntegerOutput
On exit: info=0 unless the routine detects an error (see Section 6).

6
Error Indicators and Warnings

info<0
If info=-i, argument i had an illegal value. An explanatory message is output, and execution of the program is terminated.
info=1
The Jacobi-type procedure failed to converge.

7
Accuracy

The computed generalized singular value decomposition is nearly the exact generalized singular value decomposition for nearby matrices A+E  and B+F , where
E2 = Oε A2 ​ and ​ F2 = Oε B2 ,  
and ε  is the machine precision. See Section 4.12 of Anderson et al. (1999) for further details.

8
Parallelism and Performance

f08vaf (dggsvd) makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the X06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this routine. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9
Further Comments

The complex analogue of this routine is f08vnf (zggsvd).

10
Example

This example finds the generalized singular value decomposition
A = U Σ1 0R QT ,   B = V Σ2 0R QT ,  
where
A = 1 2 3 3 2 1 4 5 6 7 8 8   and   B = -2 -3 3 4 6 5 ,  
together with estimates for the condition number of R and the error bound for the computed generalized singular values.
The example program assumes that mn, and would need slight modification if this is not the case.

10.1
Program Text

Program Text (f08vafe.f90)

10.2
Program Data

Program Data (f08vafe.d)

10.3
Program Results

Program Results (f08vafe.r)