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_dgejsv (f08kh)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_lapack_dgejsv (f08kh) computes the singular value decomposition (SVD) of a real m by n matrix A where mn, and optionally computes the left and/or right singular vectors. nag_lapack_dgejsv (f08kh) implements the preconditioned Jacobi SVD of Drmac and Veselic. This is the expert driver function that calls nag_lapack_dgesvj (f08kj) after certain preconditioning. In most cases nag_lapack_dgesvd (f08kb) or nag_lapack_dgesdd (f08kd) is sufficient to obtain the SVD of a real matrix. These are much simpler to use and also handle the case m<n.

Syntax

[a, sva, u, v, work, iwork, info] = f08kh(joba, jobu, jobv, jobr, jobt, jobp, a, 'm', m, 'n', n)
[a, sva, u, v, work, iwork, info] = nag_lapack_dgejsv(joba, jobu, jobv, jobr, jobt, jobp, a, 'm', m, 'n', n)

Description

The SVD is written as
A = UΣVT ,  
where Σ is an m by n matrix which is zero except for its n diagonal elements, U is an m by m orthogonal matrix, and V is an n by n orthogonal matrix. The diagonal elements of Σ are the singular values of A in descending order of magnitude. The columns of U and V are the left and the right singular vectors of A. The diagonal of Σ is computed and stored in the array sva.

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
Drmac Z and Veselic K (2008a) New fast and accurate Jacobi SVD algorithm I SIAM J. Matrix Anal. Appl. 29 4
Drmac Z and Veselic K (2008b) New fast and accurate Jacobi SVD algorithm II SIAM J. Matrix Anal. Appl. 29 4
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore

Parameters

Compulsory Input Parameters

1:     joba – string (length ≥ 1)
Specifies the form of pivoting for the QR factorization stage; whether an estimate of the condition number of the scaled matrix is required; and the form of rank reduction that is performed.
joba='C'
The initial QR factorization of the input matrix is performed with column pivoting; no estimate of condition number is computed; and, the rank is reduced by only the underflowed part of the triangular factor R. This option works well (high relative accuracy) if A=BD, with well-conditioned B and arbitrary diagonal matrix D. The accuracy cannot be spoiled by column scaling. The accuracy of the computed output depends on the condition of B, and the procedure aims at the best theoretical accuracy.
joba='E'
Computation as with joba='C' with an additional estimate of the condition number of B. It provides a realistic error bound.
joba='F'
The initial QR factorization of the input matrix is performed with full row and column pivoting; no estimate of condition number is computed; and, the rank is reduced by only the underflowed part of the triangular factor R. If A=D1×C×D2 with ill-conditioned diagonal scalings D1, D2, and well-conditioned matrix C, this option gives higher accuracy than the joba='C' option. If the structure of the input matrix is not known, and relative accuracy is desirable, then this option is advisable.
joba='G'
Computation as with joba='F' with an additional estimate of the condition number of B, where A=DB (i.e., B=C×D2). If A has heavily weighted rows, then using this condition number gives too pessimistic an error bound.
joba='A'
Computation as with joba='C' except in the treatment of rank reduction. In this case, small singular values are to be considered as noise and, if found, the matrix is treated as numerically rank deficient. The computed SVD A=UΣVT restores A up to fm,n×ε×A, where ε is machine precision. This gives the procedure licence to discard (set to zero) all singular values below n×ε×A.
joba='R'
Similar to joba='A'. The rank revealing property of the initial QR factorization is used to reveal (using the upper triangular factor) a gap σr+1<εσr in which case the numerical rank is declared to be r. The SVD is computed with absolute error bounds, but more accurately than with joba='A'.
Constraint: joba='C', 'E', 'F', 'G', 'A' or 'R'.
2:     jobu – string (length ≥ 1)
Specifies options for computing the left singular vectors U.
jobu='U'
The first n left singular vectors (columns of U) are computed and returned in the array u.
jobu='F'
All m left singular vectors are computed and returned in the array u.
jobu='W'
No left singular vectors are computed, but the array u (with ldum and second dimension at least n) is available as workspace for computing right singular values. See the description of u.
jobu='N'
No left singular vectors are computed. u is not referenced.
Constraint: jobu='U', 'F', 'W' or 'N'.
3:     jobv – string (length ≥ 1)
Specifies options for computing the right singular vectors V.
jobv='V'
the n right singular vectors (columns of V) are computed and returned in the array v; Jacobi rotations are not explicitly accumulated.
jobv='J'
the n right singular vectors (columns of V) are computed and returned in the array v, but they are computed as the product of Jacobi rotations. This option is allowed only if jobu='U' or 'F', i.e., in computing the full SVD.
jobv='W'
No right singular values are computed, but the array v (with ldvn and second dimension at least n) is available as workspace for computing left singular values. See the description of v.
jobv='N'
No right singular vectors are computed. v is not referenced.
Constraints:
  • jobv='V', 'J', 'W' or 'N';
  • if jobu='W' or 'N', jobv'J'.
4:     jobr – string (length ≥ 1)
Suggested value: jobr='R'.
Specifies the conditions under which columns of A are to be set to zero. This effectively specifies a lower limit on the range of singular values; any singular values below this limit are (through column zeroing) set to zero. If A0 is scaled so that the largest column (in the Euclidean norm) of cA is equal to the square root of the overflow threshold, then jobr allows the function to kill columns of A whose norm in cA is less than sfmin (for jobr='R'), or less than sfmin/ε (otherwise). sfmin is the safe range argument, as returned by function nag_machine_real_safe (x02am).
jobr='N'
Only set to zero those columns of A for which the norm of corresponding column of cA<sfmin/ε, that is, those columns that are effectively zero (to machine precision) anyway. If the condition number of A is greater than the overflow threshold λ, where λ is the value returned by nag_machine_real_largest (x02al), you are recommended to use function nag_lapack_dgesvj (f08kj).
jobr='R'
Set to zero those columns of A for which the norm of the corresponding column of cA<sfmin. This approximately represents a restricted range for σcA of sfmin,λ.
For computing the singular values in the full range from the safe minimum up to the overflow threshold use nag_lapack_dgesvj (f08kj)
Constraint: jobr='N' or 'R'.
5:     jobt – string (length ≥ 1)
Specifies, in the case n=m, whether the function is permitted to use the transpose of A for improved efficiency. If the matrix is square then the procedure may use transposed A if AT seems to be better with respect to convergence. If the matrix is not square, jobt is ignored. The decision is based on two values of entropy over the adjoint orbit of ATA. See the descriptions of work6 and work7.
jobt='T'
If n=m, perform an entropy test and then transpose if the test indicates possibly faster convergence of the Jacobi process if AT is taken as input. If A is replaced with AT, then the row pivoting is included automatically.
jobt='N'
No entropy test and no transposition is performed.
The option jobt='T' can be used to compute only the singular values, or the full SVD (U, Σ and V). In the case where only one set of singular vectors (U or V) is required, the caller must still provide both u and v, as one of the matrices is used as workspace if the matrix A is transposed. See the descriptions of u and v
Constraint: jobt='T' or 'N'.
6:     jobp – string (length ≥ 1)
Specifies whether the function should be allowed to introduce structured perturbations to drown denormalized numbers. For details see Drmac and Veselic (2008a) and Drmac and Veselic (2008b). For the sake of simplicity, these perturbations are included only when the full SVD or only the singular values are requested.
jobp='P'
Introduce perturbation if A is found to be very badly scaled (introducing denormalized numbers).
jobp='N'
Do not perturb.
Constraint: jobp='P' or 'N'.
7:     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.

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:     n int64int32nag_int scalar
Default: the second dimension of the array a.
n, the number of columns of the matrix A.
Constraint: mn0.

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.
The contents of a are overwritten.
2:     svan – double array
The, possibly scaled, singular values of A.
The singular values of A are σi=αsvai, for i=1,2,,n, where α=work1/work2. Normally α=1 and no scaling is required to obtain the singular values. However, if the largest singular value of A overflows or if small singular values have been saved from underflow by scaling the input matrix A, then α1.
If jobr='R' then some of the singular values may be returned as exact zeros because they are below the numerical rank threshold or are denormalized numbers.
3:     uldu: – double array
The first dimension, ldu, of the array u will be
  • if jobu='F', 'U' or 'W', ldu=max1,m;
  • otherwise ldu=1.
The second dimension of the array u will be max1,m if jobu='F', max1,n if jobu='U' or 'W' and 1 otherwise.
If jobu='U', u contains the m by n matrix of the left singular vectors.
If jobu='F', u contains the m by m matrix of the left singular vectors, including an orthonormal basis of the orthogonal complement of Range(A).
If jobu='W' and (jobv='V' and jobt='T' and m=n), then u is used as workspace if the procedure replaces A with AT. In that case, V is computed in u as left singular vectors of AT and then copied back to the array v.
If jobu='N', u is not referenced.
4:     vldv: – double array
The first dimension, ldv, of the array v will be
  • if jobv='V', 'J' or 'W', ldv=max1,n;
  • otherwise ldv=1.
The second dimension of the array v will be max1,n if jobv='V', 'J' or 'W' and 1 otherwise.
If jobv='V' or 'J', v contains the n by n matrix of the right singular vectors.
If jobv='W' and (jobu='U' and jobt='T' and m=n), then v is used as workspace if the procedure replaces A with AT. In that case, U is computed in v as right singular vectors of AT and then copied back to the array u.
If jobv='N', v is not referenced.
5:     worklwork – double array
Contains information about the completed job.
work1
α=work1/work2 is the scaling factor such that σi=αsvai, for i=1,2,,n are the computed singular values of A. (See the description of sva.)
work2
See the description of work1.
work3
sconda, an estimate for the condition number of column equilibrated A (if joba='E' or 'G'). sconda is an estimate of RTR-11. It is computed using nag_lapack_dpocon (f07fg). It satisfies n-14×scondaR-12n14×sconda where R is the triangular factor from the QR factorization of A. However, if R is truncated and the numerical rank is determined to be strictly smaller than n, sconda is returned as -1, thus indicating that the smallest singular values might be lost.
If full SVD is needed, and you are familiar with the details of the method, the following two condition numbers are useful for the analysis of the algorithm.
work4
An estimate of the scaled condition number of the triangular factor in the first QR factorization.
work5
An estimate of the scaled condition number of the triangular factor in the second QR factorization.
The following two parameters are computed if jobt='T'.
work6
The entropy of ATA: this is the Shannon entropy of diagATA/traceATA taken as a point in the probability simplex.
work7
The entropy of AAT.
6:     iworkm+3×n int64int32nag_int array
Contains information about the completed job.
iwork1
The numerical rank of A determined after the initial QR factorization with pivoting. See the descriptions of joba and jobr.
iwork2
The number of computed nonzero singular values.
iwork3
If nonzero, a warning message: If iwork3=1 then some of the column norms of A were denormalized (tiny) numbers. The requested high accuracy is not warranted by the data.
7:     info int64int32nag_int scalar
info=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

Cases prefixed with W are classified as warnings and do not generate an error of type NAG:error_n. See nag_issue_warnings.

   info<0
If info=-i, argument i had an illegal value. An explanatory message is output, and execution of the program is terminated.
W  info>0
nag_lapack_dgejsv (f08kh) did not converge in the allowed number of iterations (30). The computed values might be inaccurate.

Accuracy

The computed singular value decomposition is nearly the exact singular value decomposition for a nearby matrix A+E , where
E2 = Oε A2 ,  
and ε  is the machine precision. In addition, the computed singular vectors are nearly orthogonal to working precision. See Section 4.9 of Anderson et al. (1999) for further details.

Further Comments

nag_lapack_dgejsv (f08kh) implements a preconditioned Jacobi SVD algorithm. It uses nag_lapack_dgeqrf (f08ae), nag_lapack_dgelqf (f08ah) and nag_lapack_dgeqp3 (f08bf) as preprocessors and preconditioners. Optionally, an additional row pivoting can be used as a preprocessor, which in some cases results in much higher accuracy. An example is matrix A with the structure A=D1CD2, where D1, D2 are arbitrarily ill-conditioned diagonal matrices and C is a well-conditioned matrix. In that case, complete pivoting in the first QR factorizations provides accuracy dependent on the condition number of C, and independent of D1, D2. Such higher accuracy is not completely understood theoretically, but it works well in practice. Further, if A can be written as A=BD, with well-conditioned B and some diagonal D, then the high accuracy is guaranteed, both theoretically and in software, independent of D.

Example

This example finds the singular values and left and right singular vectors of the 6 by 4 matrix
A = 2.27 -1.54 1.15 -1.94 0.28 -1.67 0.94 -0.78 -0.48 -3.09 0.99 -0.21 1.07 1.22 0.79 0.63 -2.35 2.93 -1.45 2.30 0.62 -7.39 1.03 -2.57 ,  
together with the condition number of A and approximate error bounds for the computed singular values and vectors.
function f08kh_example


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

% General 6-by-4 real matrix A
m = int64(6);
n = int64(4);
a = [2.27, -1.54,  1.15, -1.94;
     0.28, -1.67,  0.94, -0.78;
    -0.48, -3.09,  0.99, -0.21;
     1.07,  1.22,  0.79,  0.63;
    -2.35,  2.93, -1.45,  2.30;
     0.62, -7.39,  1.03, -2.57];

% Compute the SVD of A (A = U*S*V^T, m >= n)
joba = 'estimated condition number';
jobu = 'U required';
jobv = 'V required';
jobr = 'restricted range';
jobt = 'No entropy test';
jobp = 'No perturbations';
[~, s, u, v, work, iwork, info] = ...
  f08kh( ...
         joba, jobu, jobv, jobr, jobt, jobp, a);

% Approximate error bound for S, s(1) = norm(A).
eps = x02aj;
serrbd = eps*s(1);

% Print solution
if (abs(work(1)-work(2)) < 2*eps)
  % No scaling required
  fprintf('Singular values:\n');
  disp(transpose(s));
else
  fprintf('Scaled singular values:\n');
  disp(transpose(s));
  fprintf('\nFor true singular values, multiply by a/b\n');
  fprintf('   where a=%13.5e and b=%13.5e.\n', work(1), work(2));
end

[ifail] = x04ca( ...
                 'Gen', ' ', u, 'Left singular vectors');
fprintf('\n');
[ifail] = x04ca( ...
                 'Gen', ' ', v, 'Right singular vectors');

% Estimate reciprocal condition numbers for S
[rcondu, info] = f08fl( ...
                        'Left', m, n, s);
[rcondv, info] = f08fl( ...
                        'Right', m, n, s);

% display approximate error bounds
fprintf('\nEstimate of the condition number of column equilibrated A\n');
fprintf('%11.1e\n', work(3));
fprintf('\nError estimate  for S:\n');
fprintf('%11.1e\n', serrbd);
fprintf('\nError estimates for U:\n');
fprintf('%11.1e ',serrbd./rcondu);
fprintf('\n\nError estimates for V:\n');
fprintf('%11.1e ',serrbd./rcondv);
fprintf('\n');


f08kh example results

Singular values:
    9.9966    3.6831    1.3569    0.5000

 Left singular vectors
          1       2       3       4
 1   0.2774 -0.6003 -0.1277  0.1323
 2   0.2020 -0.0301  0.2805  0.7034
 3   0.2918  0.3348  0.6453  0.1906
 4  -0.0938 -0.3699  0.6781 -0.5399
 5  -0.4213  0.5266  0.0413 -0.0575
 6   0.7816  0.3353 -0.1645 -0.3957

 Right singular vectors
          1       2       3       4
 1   0.1921 -0.8030  0.0041 -0.5642
 2  -0.8794 -0.3926 -0.0752  0.2587
 3   0.2140 -0.2980  0.7827  0.5027
 4  -0.3795  0.3351  0.6178 -0.6017

Estimate of the condition number of column equilibrated A
    9.0e+00

Error estimate  for S:
    1.1e-15

Error estimates for U:
    1.8e-16     4.8e-16     1.3e-15     2.2e-15 

Error estimates for V:
    1.8e-16     4.8e-16     1.3e-15     1.3e-15 

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