NAG FL Interface
f08vcf (dggsvd3)

Settings help

FL Name Style:


FL Specification Language:


1 Purpose

f08vcf computes the generalized singular value decomposition (GSVD) of an m×n real matrix A and a p×n real matrix B.

2 Specification

Fortran Interface
Subroutine f08vcf ( jobu, jobv, jobq, m, n, p, k, l, a, lda, b, ldb, alpha, beta, u, ldu, v, ldv, q, ldq, work, lwork, iwork, info)
Integer, Intent (In) :: m, n, p, lda, ldb, ldu, ldv, ldq, lwork
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(1,lwork))
Character (1), Intent (In) :: jobu, jobv, jobq
C Header Interface
#include <nag.h>
void  f08vcf_ (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[], const Integer *lwork, Integer iwork[], Integer *info, const Charlen length_jobu, const Charlen length_jobv, const Charlen length_jobq)
The routine may be called by the names f08vcf, nagf_lapackeig_dggsvd3 or its LAPACK name dggsvd3.

3 Description

Given an m×n real matrix A and a p×n real matrix B, 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 l be the effective numerical rank of B and (k+l) be the effective numerical rank of the matrix ( A B ) , then the first k generalized singular values are infinite and the remaining l are finite. R is a (k+l)×(k+l) nonsingular upper triangular matrix, D1 and D2 are m×(k+l) and p×(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 (R33)ij 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×n nonsingular matrix, then the GSVD of A and B implicitly gives the SVD of AB-1:
A B-1 = U (D1D2−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 setting
X = Q ( I 0 0 R-1 ) .  
A two stage process is used to compute the GSVD of the matrix pair (A,B). The pair is first reduced to upper triangular form by orthogonal transformations using f08vgf. The GSVD of the resulting upper triangular matrix pair is then performed by f08yef which uses a variant of the Kogbetliantz algorithm (a cyclic Jacobi method).

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 https://www.netlib.org/lapack/lug
Golub G H and Van Loan C F (2012) Matrix Computations (4th 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 Integer Input
On entry: m, the number of rows of the matrix A.
Constraint: m0.
5: n Integer Input
On entry: n, the number of columns of the matrices A and B.
Constraint: n0.
6: p Integer Input
On entry: p, the number of rows of the matrix B.
Constraint: p0.
7: k Integer Output
8: l Integer Output
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: a(lda,*) Real (Kind=nag_wp) array Input/Output
Note: the second dimension of the array a must be at least max(1,n).
On entry: the m×n matrix A.
On exit: contains the triangular matrix R, or part of R. See Section 3 for details.
10: lda Integer Input
On entry: the first dimension of the array a as declared in the (sub)program from which f08vcf is called.
Constraint: ldamax(1,m).
11: b(ldb,*) Real (Kind=nag_wp) array Input/Output
Note: the second dimension of the array b must be at least max(1,n).
On entry: the p×n matrix B.
On exit: contains the triangular matrix R if m-k-l<0. See Section 3 for details.
12: ldb Integer Input
On entry: the first dimension of the array b as declared in the (sub)program from which f08vcf is called.
Constraint: ldbmax(1,p).
13: alpha(n) Real (Kind=nag_wp) array Output
On exit: see the description of beta.
14: beta(n) Real (Kind=nag_wp) array Output
On exit: alpha and beta contain the generalized singular value pairs of A and B, αi and βi ;
  • alpha(1:k) = 1 ,
  • beta(1:k) = 0 ,
and if m-k-l0 ,
  • 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 .
The notation alpha(k:n) above refers to consecutive elements alpha(i), for i=k,,n.
15: u(ldu,*) Real (Kind=nag_wp) array Output
Note: the second dimension of the array u must be at least max(1,m) if jobu='U', and at least 1 otherwise.
On exit: if jobu='U', u contains the m×m orthogonal matrix U.
If jobu='N', u is not referenced.
16: ldu Integer Input
On entry: the first dimension of the array u as declared in the (sub)program from which f08vcf is called.
Constraints:
  • if jobu='U', ldu max(1,m) ;
  • otherwise ldu1.
17: v(ldv,*) Real (Kind=nag_wp) array Output
Note: the second dimension of the array v must be at least max(1,p) if jobv='V', and at least 1 otherwise.
On exit: if jobv='V', v contains the p×p orthogonal matrix V.
If jobv='N', v is not referenced.
18: ldv Integer Input
On entry: the first dimension of the array v as declared in the (sub)program from which f08vcf is called.
Constraints:
  • if jobv='V', ldv max(1,p) ;
  • otherwise ldv1.
19: q(ldq,*) Real (Kind=nag_wp) array Output
Note: the second dimension of the array q must be at least max(1,n) if jobq='Q', and at least 1 otherwise.
On exit: if jobq='Q', q contains the n×n orthogonal matrix Q.
If jobq='N', q is not referenced.
20: ldq Integer Input
On entry: the first dimension of the array q as declared in the (sub)program from which f08vcf is called.
Constraints:
  • if jobq='Q', ldq max(1,n) ;
  • otherwise ldq1.
21: work(max(1,lwork)) Real (Kind=nag_wp) array Workspace
On exit: if info=0, work(1) contains the minimum value of lwork required for optimal performance.
22: lwork Integer Input
On entry: the dimension of the array work as declared in the (sub)program from which f08vcf is called.
If lwork=−1, 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.
Suggested value: for optimal performance, lwork must generally be larger than the minimum; increase workspace by, say, nb×(n+1), where nb is the optimal block size.
Constraints:
if lwork−1,
  • if jobv='V', lworkmax(3*n+1,m,p);
  • if jobv='N', lworkmax(3*n+1,m).
23: iwork(n) Integer array Output
On exit: stores the sorting information. More precisely, if I is the ordered set of indices of alpha containing C (denote as alpha(I), see beta), then the corresponding elements iwork(I) contain the swap pivots, J, that sorts I such that alpha(I) is in descending numerical order.
The following pseudocode sorts the set I:
for ​iI j=Ji swap ​Ii​ and ​Ij end
24: info Integer Output
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

f08vcf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08vcf 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

This routine replaces the deprecated routine f08vaf which used an unblocked algorithm and, therefore, did not make best use of Level 3 BLAS routines.
The complex analogue of this routine is f08vqf.

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 (f08vcfe.f90)

10.2 Program Data

Program Data (f08vcfe.d)

10.3 Program Results

Program Results (f08vcfe.r)