NAG Library Routine Document

f08vqf  (zggsvd3)

 Contents

    1  Purpose
    7  Accuracy

1
Purpose

f08vqf (zggsvd3) computes the generalized singular value decomposition (GSVD) of an m by n complex matrix A and a p by n complex matrix B.

2
Specification

Fortran Interface
Subroutine f08vqf ( jobu, jobv, jobq, m, n, p, k, l, a, lda, b, ldb, alpha, beta, u, ldu, v, ldv, q, ldq, work, lwork, rwork, 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 (Out):: alpha(n), beta(n), rwork(2*n)
Complex (Kind=nag_wp), Intent (Inout):: a(lda,*), b(ldb,*), u(ldu,*), v(ldv,*), q(ldq,*)
Complex (Kind=nag_wp), Intent (Out):: work(max(1,lwork))
Character (1), Intent (In):: jobu, jobv, jobq
C Header Interface
#include nagmk26.h
void  f08vqf_ ( const char *jobu, const char *jobv, const char *jobq, const Integer *m, const Integer *n, const Integer *p, Integer *k, Integer *l, Complex a[], const Integer *lda, Complex b[], const Integer *ldb, double alpha[], double beta[], Complex u[], const Integer *ldu, Complex v[], const Integer *ldv, Complex q[], const Integer *ldq, Complex work[], const Integer *lwork, double rwork[], 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 zggsvd3.

3
Description

Given an m by n complex matrix A and a p by n complex matrix B, the generalized singular value decomposition is given by
UH A Q = D1 0 R ,   VH B Q = D2 0 R ,  
where U, V and Q are unitary 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 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 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×B-1:
A B-1 = U D1 D2-1 VH .  
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:
AH Ax=λ BH Bx .  
In some literature, the GSVD of A and B is presented in the form
UH A X = 0D1 ,   VH 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 unitary transformations using f08vuf (zggsvp3). The GSVD of the resulting upper triangular matrix pair is then performed by f08ysf (ztgsja) 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 http://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 unitary 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 unitary 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 unitary 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* – Complex (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 f08vqf (zggsvd3) is called.
Constraint: ldamax1,m.
11:   bldb* – Complex (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 f08vqf (zggsvd3) 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* – Complex (Kind=nag_wp) arrayOutput
Note: the second dimension of the array u must be at least max1,m if jobu='U', and at least 1 otherwise.
On exit: if jobu='U', u contains the m by m unitary 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 f08vqf (zggsvd3) is called.
Constraints:
  • if jobu='U', ldu max1,m ;
  • otherwise ldu1.
17:   vldv* – Complex (Kind=nag_wp) arrayOutput
Note: the second dimension of the array v must be at least max1,p if jobv='V', and at least 1 otherwise.
On exit: if jobv='V', v contains the p by p unitary 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 f08vqf (zggsvd3) is called.
Constraints:
  • if jobv='V', ldv max1,p ;
  • otherwise ldv1.
19:   qldq* – Complex (Kind=nag_wp) arrayOutput
Note: the second dimension of the array q must be at least max1,n if jobq='Q', and at least 1 otherwise.
On exit: if jobq='Q', q contains the n by n unitary 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 f08vqf (zggsvd3) is called.
Constraints:
  • if jobq='Q', ldq max1,n ;
  • otherwise ldq1.
21:   workmax1,lwork – Complex (Kind=nag_wp) arrayWorkspace
On exit: if info=0, the real part of work1 contains the minimum value of lwork required for optimal performance.
22:   lwork – IntegerInput
On entry: the dimension of the array work as declared in the (sub)program from which f08vqf (zggsvd3) 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 jobv='V', lworkmaxn+1,m,p;
  • if jobv='N', lworkmaxn+1,m.
23:   rwork2×n – Real (Kind=nag_wp) arrayWorkspace
24:   iworkn – Integer arrayOutput
On exit: stores the sorting information. More precisely, if I is the ordered set of indices of alpha containing C (denote as alphaI, see beta), then the corresponding elements iworkI contain the swap pivots, J, that sorts I such that alphaI is in descending numerical order.
The following pseudocode sorts the set I:
for ​iI j=Ji swap ​Ii​ and ​Ij end  
25:   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

f08vqf (zggsvd3) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08vqf (zggsvd3) 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 f08vnf (zggsvd) which used an unblocked algorithm and therefore did not make best use of level 3 BLAS routines.
The diagonal elements of the matrix R are real.
The real analogue of this routine is f08vcf (dggsvd3).

10
Example

This example finds the generalized singular value decomposition
A = U Σ1 0R QH ,   B = V Σ2 0R QH ,  
where
A = 0.96-0.81i -0.03+0.96i -0.91+2.06i -0.05+0.41i -0.98+1.98i -1.20+0.19i -0.66+0.42i -0.81+0.56i 0.62-0.46i 1.01+0.02i 0.63-0.17i -1.11+0.60i 0.37+0.38i 0.19-0.54i -0.98-0.36i 0.22-0.20i 0.83+0.51i 0.20+0.01i -0.17-0.46i 1.47+1.59i 1.08-0.28i 0.20-0.12i -0.07+1.23i 0.26+0.26i  
and
B = 1 0 -1 0 0 1 0 -1 ,  
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 (f08vqfe.f90)

10.2
Program Data

Program Data (f08vqfe.d)

10.3
Program Results

Program Results (f08vqfe.r)

© The Numerical Algorithms Group Ltd, Oxford, UK. 2017