PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_lapack_dggsvp (f08ve)
Purpose
nag_lapack_dggsvp (f08ve) uses orthogonal transformations to simultaneously reduce the by matrix and the by matrix 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
,
and
such that
where the
by
matrix
and
by
matrix
are nonsingular upper triangular;
is
by
upper triangular if
and is
by
upper trapezoidal otherwise.
is the effective numerical rank of the
by
matrix
.
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:
– string (length ≥ 1)
-
If
, the orthogonal matrix
is computed.
If , is not computed.
Constraint:
or .
- 2:
– string (length ≥ 1)
-
If
, the orthogonal matrix
is computed.
If , is not computed.
Constraint:
or .
- 3:
– string (length ≥ 1)
-
If
, the orthogonal matrix
is computed.
If , is not computed.
Constraint:
or .
- 4:
– double array
-
The first dimension of the array
a must be at least
.
The second dimension of the array
a must be at least
.
The by matrix .
- 5:
– double array
-
The first dimension of the array
b must be at least
.
The second dimension of the array
b must be at least
.
The by matrix .
- 6:
– double scalar
- 7:
– double scalar
-
tola and
tolb are the thresholds to determine the effective numerical rank of matrix
and a subblock of
. Generally, they are set to
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:
– int64int32nag_int scalar
-
Default:
the first dimension of the array
a.
, the number of rows of the matrix .
Constraint:
.
- 2:
– int64int32nag_int scalar
-
Default:
the first dimension of the array
b.
, the number of rows of the matrix .
Constraint:
.
- 3:
– int64int32nag_int scalar
-
Default:
the second dimension of the array
a.
, the number of columns of the matrices and .
Constraint:
.
Output Parameters
- 1:
– double array
-
The first dimension of the array
a will be
.
The second dimension of the array
a will be
.
Contains the triangular (or trapezoidal) matrix described in
Description.
- 2:
– double array
-
The first dimension of the array
b will be
.
The second dimension of the array
b will be
.
Contains the triangular matrix described in
Description.
- 3:
– int64int32nag_int scalar
- 4:
– int64int32nag_int scalar
-
k and
l specify the dimension of the subblocks
and
as described in
Description;
is the effective numerical rank of
.
- 5:
– double array
-
The first dimension,
, of the array
u will be
- if , ;
- otherwise .
The second dimension of the array
u will be
if
and
otherwise.
If
,
u contains the orthogonal matrix
.
If
,
u is not referenced.
- 6:
– double array
-
The first dimension,
, of the array
v will be
- if , ;
- otherwise .
The second dimension of the array
v will be
if
and
otherwise.
If
,
v contains the orthogonal matrix
.
If
,
v is not referenced.
- 7:
– double array
-
The first dimension,
, of the array
q will be
- if , ;
- otherwise .
The second dimension of the array
q will be
if
and
otherwise.
If
,
q contains the orthogonal matrix
.
If
,
q is not referenced.
- 8:
– int64int32nag_int scalar
unless the function detects an error (see
Error Indicators and Warnings).
Error Indicators and Warnings
-
If , parameter 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
and
, where
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
of the matrix pair
, where
Open in the MATLAB editor:
f08ve_example
function f08ve_example
fprintf('f08ve example results\n\n');
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];
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);
[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)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015