hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_lapack_dggsvp (f08ve)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_lapack_dggsvp (f08ve) uses orthogonal transformations to simultaneously reduce the m by n matrix A and the p by n matrix B to upper triangular form. This factorization is usually used as a preprocessing step for computing the generalized singular value decomposition (GSVD).

Syntax

[a, b, k, l, u, v, q, info] = f08ve(jobu, jobv, jobq, a, b, tola, tolb, 'm', m, 'p', p, 'n', n)
[a, b, k, l, u, v, q, info] = nag_lapack_dggsvp(jobu, jobv, jobq, a, b, tola, tolb, 'm', m, 'p', p, 'n', n)

Description

nag_lapack_dggsvp (f08ve) computes orthogonal matrices U, V and Q such that
UTAQ= n-k-lklk0A12A13l00A23m-k-l000() , if ​m-k-l0; n-k-lklk0A12A13m-k00A23() , if ​m-k-l<0; VTBQ= n-k-lkll00B13p-l000()  
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-l0 and is m-k by l upper trapezoidal otherwise. k+l is the effective numerical rank of the m+p by n matrix ATBTT.
This decomposition is usually used as the preprocessing step for computing the Generalized Singular Value Decomposition (GSVD), see function nag_lapack_dggsvd (f08va).

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

Parameters

Compulsory Input Parameters

1:     jobu – string (length ≥ 1)
If jobu='U', the orthogonal matrix U is computed.
If jobu='N', U is not computed.
Constraint: jobu='U' or 'N'.
2:     jobv – string (length ≥ 1)
If jobv='V', the orthogonal matrix V is computed.
If jobv='N', V is not computed.
Constraint: jobv='V' or 'N'.
3:     jobq – string (length ≥ 1)
If jobq='Q', the orthogonal matrix Q is computed.
If jobq='N', Q is not computed.
Constraint: jobq='Q' or 'N'.
4:     alda: – double array
The first dimension of the array a must be at least max1,m.
The second dimension of the array a must be at least max1,n.
The m by n matrix A.
5:     bldb: – double array
The first dimension of the array b must be at least max1,p.
The second dimension of the array b must be at least max1,n.
The p by n matrix B.
6:     tola – double scalar
7:     tolb – double scalar
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=maxm,nAε, tolb=maxp,nBε,  
where ε  is the machine precision.
The size of tola and tolb may affect the size of backward errors of the decomposition.

Optional Input Parameters

1:     m int64int32nag_int scalar
Default: the first dimension of the array a.
m, the number of rows of the matrix A.
Constraint: m0.
2:     p int64int32nag_int scalar
Default: the first dimension of the array b.
p, the number of rows of the matrix B.
Constraint: p0.
3:     n int64int32nag_int scalar
Default: the second dimension of the array a.
n, the number of columns of the matrices A and B.
Constraint: n0.

Output Parameters

1:     alda: – double array
The first dimension of the array a will be max1,m.
The second dimension of the array a will be max1,n.
Contains the triangular (or trapezoidal) matrix described in Description.
2:     bldb: – double array
The first dimension of the array b will be max1,p.
The second dimension of the array b will be max1,n.
Contains the triangular matrix described in Description.
3:     k int64int32nag_int scalar
4:     l int64int32nag_int scalar
k and l specify the dimension of the subblocks k and l as described in Description; k+l is the effective numerical rank of aTbTT.
5:     uldu: – double array
The first dimension, ldu, of the array u will be
  • if jobu='U', ldu= max1,m ;
  • otherwise ldu=1.
The second dimension of the array u will be max1,m if jobu='U' and 1 otherwise.
If jobu='U', u contains the orthogonal matrix U.
If jobu='N', u is not referenced.
6:     vldv: – double array
The first dimension, ldv, of the array v will be
  • if jobv='V', ldv= max1,p ;
  • otherwise ldv=1.
The second dimension of the array v will be max1,p if jobv='V' and 1 otherwise.
If jobv='V', v contains the orthogonal matrix V.
If jobv='N', v is not referenced.
7:     qldq: – double array
The first dimension, ldq, of the array q will be
  • if jobq='Q', ldq= max1,n ;
  • otherwise ldq=1.
The second dimension of the array q will be max1,n if jobq='Q' and 1 otherwise.
If jobq='Q', q contains the orthogonal matrix Q.
If jobq='N', q is not referenced.
8:     info int64int32nag_int scalar
info=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

   info=-i
If info=-i, parameter i had an illegal value on entry. The parameters are numbered as follows:
1: jobu, 2: jobv, 3: jobq, 4: m, 5: p, 6: n, 7: a, 8: lda, 9: b, 10: ldb, 11: tola, 12: tolb, 13: k, 14: l, 15: u, 16: ldu, 17: v, 18: ldv, 19: q, 20: ldq, 21: iwork, 22: tau, 23: work, 24: info.
It is possible that info refers to a parameter that is omitted from the MATLAB interface. This usually indicates that an error in one of the other input parameters has caused an incorrect value to be inferred.

Accuracy

The computed factorization is nearly the exact factorization for nearby matrices A+E and B+F, where
E2 = OεA2   and   F2= OεB2,  
and ε is the machine precision.

Further Comments

The complex analogue of this function is nag_lapack_zggsvp (f08vs).

Example

This example finds the generalized factorization
A = UΣ1 0 S QT ,   B= VΣ2 0 T QT ,  
of the matrix pair AB, where
A = 123 321 456 788   and   B= -2-33 465 .  
function f08ve_example


fprintf('f08ve example results\n\n');

% Generalized singular values of (A,B) where:
m = 4;
n = 3;
a = [ 1  2  3;
      3  2  1;
      4  5  6;
      7  8  8]; 
p = 2;
b = [-2 -3  3;
      4  6  5]; 

% Reduce A and B to upper triangular form S = U^T A Q, T = V^T B Q
tola = max(m,n)*norm(a,1)*x02aj;
tolb = max(p,n)*norm(a,1)*x02aj;
[S, T, k, l, U, V, Q, info] = ...
  f08ve( ...
	 'U', 'V', 'Q', a, b, tola, tolb);

% Compute singular values
[S, T, alpha, beta, U, V, Q, ncycle, info] = ...
  f08ye( ...
	 'U', 'V', 'Q', k, l, S, T, tola, tolb, U, V, Q);

fprintf('Number of infinite generalized singular values = %3d\n',k);
fprintf('Number of   finite generalized singular values = %3d\n',l);
fprintf('Effective rank of the matrix pair  (A^T B^T)^T = %3d\n',k+l);
fprintf('Number of cycles of the Kogbetliantz method    = %3d\n\n',ncycle);
disp('Finite generalized singular values');
disp(alpha(k+1:k+l)./beta(k+1:k+l));


f08ve example results

Number of infinite generalized singular values =   1
Number of   finite generalized singular values =   2
Effective rank of the matrix pair  (A^T B^T)^T =   3
Number of cycles of the Kogbetliantz method    =   2

Finite generalized singular values
    1.3151
    0.0802


PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015