NAG FL Interface
f10caf (svd_​rowext_​real)

Settings help

FL Name Style:


FL Specification Language:


1 Purpose

f10caf computes the singular value decomposition (SVD) of a real m×n matrix A, optionally computing the left and/or right singular vectors by using a randomised numerical linear algebra (RNLA) method.

2 Specification

Fortran Interface
Subroutine f10caf ( jobu, jobvt, m, n, a, lda, k, rtol_abs, rtol_rel, state, s, u, ldu, vt, ldvt, r, ifail)
Integer, Intent (In) :: m, n, lda, k, ldu, ldvt
Integer, Intent (Inout) :: state(*), ifail
Integer, Intent (Out) :: r
Real (Kind=nag_wp), Intent (In) :: a(lda,*), rtol_abs, rtol_rel
Real (Kind=nag_wp), Intent (Inout) :: s(k), u(ldu,*), vt(ldvt,*)
Character (1), Intent (In) :: jobu, jobvt
C Header Interface
#include <nag.h>
void  f10caf_ (const char *jobu, const char *jobvt, const Integer *m, const Integer *n, const double a[], const Integer *lda, const Integer *k, const double *rtol_abs, const double *rtol_rel, Integer state[], double s[], double u[], const Integer *ldu, double vt[], const Integer *ldvt, Integer *r, Integer *ifail, const Charlen length_jobu, const Charlen length_jobvt)
The routine may be called by the names f10caf or nagf_rnla_svd_rowext_real.

3 Description

The SVD is written as
A = UΣVT ,  
where Σ is an m×n matrix which is zero except for its min(m,n) diagonal elements, U is an m×m orthogonal matrix, and V is an n×n orthogonal matrix. The diagonal elements of Σ are the singular values of A; they are real and non-negative, and are returned in descending order. The first min(m,n) columns of U and V are the left and right singular vectors of A.
Note that the routine returns VT, not V.
If the rank of A is r<min(m,n), then σ has r nonzero elements, and only r columns of U and V are well-defined. In this case we can reduce σ to an r×r matrix, U to an m×r matrix and VT to an r×n matrix.
f10caf is designed for efficiently computing the SVD in the case r<min(m,n). The input argument k should be greater than r by a small oversampling parameter, δ, such that k=r+δ. A reasonable value for δ, to compute the SVD to within machine precision, is δ=5. The value of δ should not vary based on m or n. If r is not known then the routine can be used iteratively to refine the estimate and accuracy of the computed SVD; using a larger value of δ than necessary increases the computational cost of the routine.
As a by-product of computing the SVD, the routine estimates r.
If the input argument k is less than r the accuracy depends on the (k+1)th singular value, σk+1. See Section 7 for more details.
A call to f10caf consists of the following:
  1. 1.A random projection is applied, Y=AΩ, where Ω is an n×k matrix. (Note that the product AΩ is computed using a Fast Fourier Transform, so can be computed in O(mnlog(k)) time.) See f10daf for more details on the random projection.
  2. 2.A pivoted QR decomposition of Y is calculated (see f08bff for more details). The rank estimate is then r such that, on the diagonal of R,
    | R rr | > εa + εr × R 11 ,  
    where εa and εr are the absolute and relative error tolerances, respectively, and r is the largest diagonal index for which the above relation holds.
  3. 3.Obtain the SVD from the QR decomposition of Y (or, depending on the rank, an approximation to the SVD) of A. This is referred to as row extraction.
Further details of the randomized SVD procedure can be found in Sections 4 and 5 of Halko et al. (2011).

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 (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore
Halko N (2012) Randomized methods for computing low-rank approximations of matrices PhD thesis
Halko N, Martinsson P G and Tropp J A (2011) Finding structure with randomness: Probabilistic algorithms for constructing approximate matrix decompositions SIAM Rev. 53(2) 217–288 https://epubs.siam.org/doi/abs/10.1137/090771806

5 Arguments

1: jobu Character(1) Input
On entry: specifies options for computing part of or none of the matrix U.
jobu='S'
The first k columns of U (the left singular vectors) are returned in the array u.
jobu='N'
No columns of U (no left singular vectors) are computed.
Constraint: jobu='S' or 'N'.
2: jobvt Character(1) Input
On entry: specifies options for computing part of or none of the matrix VT.
jobvt='S'
The first k rows of VT (the right singular vectors) are returned in the array vt.
jobvt='N'
No rows of VT (no right singular vectors) are computed.
Constraint: jobvt='S' or 'N'.
3: m Integer Input
On entry: m, the number of rows of the matrix A.
Constraint: m>0.
4: n Integer Input
On entry: n, the number of columns of the matrix A.
Constraint: n>0.
5: a(lda,*) Real (Kind=nag_wp) array Input
Note: the second dimension of the array a must be at least n.
On entry: the m×n matrix A.
6: lda Integer Input
On entry: the first dimension of the array a as declared in the (sub)program from which f10caf is called.
Constraint: ldam.
7: k Integer Input
On entry: k, number of columns in random projection, Y=AΩ.
Constraint: 0<kmin(m,n).
8: rtol_abs Real (Kind=nag_wp) Input
On entry: the absolute tolerance, εa used in defining the threshold on estimating the rank of A. If rtol_abs<0.0 then 0.0 is used unless rtol_rel0.0 in which case machine precision is used.
9: rtol_rel Real (Kind=nag_wp) Input
On entry: the relative tolerance, εr used in defining the threshold on estimating the rank of A. If rtol_rel<0.0 then 0.0 is used unless rtol_abs0.0 in which case machine precision is used.
10: state(*) Integer array Communication Array
Note: the actual argument supplied must be the array state supplied to the initialization routines g05kff or g05kgf.
On entry: contains information on the selected base generator and its current state.
On exit: contains updated information on the state of the generator.
11: s(k) Real (Kind=nag_wp) array Output
On exit: the first r elements of s contain the r largest singular values of A in descending order. The remaining values are set to zero.
12: u(ldu,*) Real (Kind=nag_wp) array Output
Note: the second dimension of the array u must be at least k if jobu='S'.
On exit: if jobu='S', u contains the first r columns of U (the left singular vectors, stored column-wise); the remaining elements of u are set to zero.
If jobu='N', u is not referenced.
13: ldu Integer Input
On entry: the first dimension of the array u as declared in the (sub)program from which f10caf is called.
Constraint: if jobu='S', ldum.
14: vt(ldvt,*) Real (Kind=nag_wp) array Output
Note: the second dimension of the array vt must be at least n if jobvt='S', and at least 1 otherwise.
On exit: if jobvt='S', vt contains the first r rows of VT (the right singular vectors); the remaining elements of vt are set to zero.
If jobvt='N', vt is not referenced.
15: ldvt Integer Input
On entry: the first dimension of the array vt as declared in the (sub)program from which f10caf is called.
Constraint: if jobvt='S', ldvtk.
16: r Integer Output
On exit: r, contains estimated rank of array a.
17: ifail Integer Input/Output
On entry: ifail must be set to 0, -1 or 1 to set behaviour on detection of an error; these values have no effect when no error is detected.
A value of 0 causes the printing of an error message and program execution will be halted; otherwise program execution continues. A value of -1 means that an error message is printed while a value of 1 means that it is not.
If halting is not appropriate, the value -1 or 1 is recommended. If message printing is undesirable, then the value 1 is recommended. Otherwise, the value 0 is recommended. When the value -1 or 1 is used it is essential to test the value of ifail on exit.
On exit: ifail=0 unless the routine detects an error or a warning has been flagged (see Section 6).

6 Error Indicators and Warnings

If on entry ifail=0 or -1, explanatory error messages are output on the current error message unit (as defined by x04aaf).
Errors or warnings detected by the routine:
ifail=1
On entry, jobu=value.
Constraint: jobu='S' or 'N'.
ifail=2
On entry, jobvt=value.
Constraint: jobvt='S' or 'N'.
ifail=3
On entry, m=value.
Constraint: m>0.
ifail=4
On entry, n=value.
Constraint: n>0.
ifail=6
On entry, lda=value.
Constraint: ldam.
ifail=7
On entry, k=value.
Constraint: 0<kmin(m,n).
ifail=8
On entry, state vector has been corrupted or not initialized.
ifail=9
On entry, jobu=value and ldu=value.
Constraint: if jobu='S', ldum.
ifail=10
On entry, jobvt=value and ldvt=value.
Constraint: if jobvt='S', ldvtk.
ifail=21
On exit, r=k, the rank of A may be larger than r.
Increase k to obtain a more accurate rank estimate.
Smallest diagonal element of R, from QR of Y, =value.
Tolerance used to determine rank =value.
ifail=22
A has effective rank of zero.
First diagonal element of R, from QR of Y, =value.
Tolerance used to determine rank =value.
ifail=-99
An unexpected error has been triggered by this routine. Please contact NAG.
See Section 7 in the Introduction to the NAG Library FL Interface for further information.
ifail=-399
Your licence key may have expired or may not have been installed correctly.
See Section 8 in the Introduction to the NAG Library FL Interface for further information.
ifail=-999
Dynamic memory allocation failed.
See Section 9 in the Introduction to the NAG Library FL Interface for further information.

7 Accuracy

The error is approximately,
A-UσV* < [1+ 1 + 4k (n-k) ] ε ,  
where,
ε= 1+7n/k · σk+1 .  
The norm on the left-hand side of the first equation is the spectral norm, and σk+1 is the (k+1)th singular value of A. More details on the error bound can be found in Sections 5 and 11 of Halko et al. (2011).

8 Parallelism and Performance

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

The total number of floating-point operations is O(log(k)mn+k2(m+n)). The first term corresponds to applying the random projection, i.e., computing Y=AΩ. The second term corresponds to the QR decomposition of Y and the steps required to obtain the SVD of the original matrix A.
Deterministic SVD solvers, such as f08kbf, require O(mn2) operations when mn and O(nm2) operations when nm.
The default values for rtol_abs and rtol_rel assume that you need an accurate approximation to A. If you only need to use a small number of singular values or singular vectors, larger values for these tolerances are appropriate. Increasing tolerances sufficiently will decrease r, the estimated rank. Decreasing r means that k can then be decreased to reduce the run-time of the routine.

10 Example

This example finds the singular values, the left and right singular vectors, and the rank of the 6×6 matrix
A = ( 2.27 0.28 -0.48 1.07 -2.35 0.62 -1.54 -1.67 -3.09 1.22 2.93 -7.39 1.15 0.94 0.99 0.79 -1.45 1.03 -1.94 -0.78 -0.21 0.63 2.30 -2.57 -1.94 -0.78 -0.21 0.63 2.30 -2.57 -1.94 -0.78 -0.21 0.63 2.30 -2.57 ) ,  
using the randomised solver, f10caf, and a deterministic solver, f08kbf for comparison.

10.1 Program Text

Program Text (f10cafe.f90)

10.2 Program Data

Program Data (f10cafe.d)

10.3 Program Results

Program Results (f10cafe.r)