f10caf computes the singular value decomposition (SVD) of a real matrix , optionally computing the left and/or right singular vectors by using a randomised numerical linear algebra (RNLA) method.
The routine may be called by the names f10caf or nagf_rnla_svd_rowext_real.
3Description
The SVD is written as
where is an matrix which is zero except for its diagonal elements, is an orthogonal matrix, and is an orthogonal matrix. The diagonal elements of are the singular values of ; they are real and non-negative, and are returned in descending order. The first columns of and are the left and right singular vectors of .
Note that the routine returns , not .
If the rank of is , then has nonzero elements, and only columns of and are well-defined. In this case we can reduce to an matrix, to an matrix and to an matrix.
f10caf is designed for efficiently computing the SVD in the case . The input argument should be greater than by a small oversampling parameter, , such that . A reasonable value for , to compute the SVD to within machine precision, is . The value of should not vary based on or . If 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 .
If the input argument is less than the accuracy depends on the th singular value, . See Section 7 for more details.
A call to f10caf consists of the following:
1.A random projection is applied, , where is an matrix. (Note that the product is computed using a Fast Fourier Transform, so can be computed in time.) See f10daf for more details on the random projection.
2.A pivoted decomposition of is calculated (see f08bff for more details). The rank estimate is then such that, on the diagonal of ,
where and are the absolute and relative error tolerances, respectively, and is the largest diagonal index for which the above relation holds.
3.Obtain the SVD from the decomposition of (or, depending on the rank, an approximation to the SVD) of . 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).
4References
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
5Arguments
1: – Character(1)Input
On entry: specifies options for computing part of or none of the matrix .
The first columns of (the left singular vectors) are returned in the array u.
No columns of (no left singular vectors) are computed.
Constraint:
or .
2: – Character(1)Input
On entry: specifies options for computing part of or none of the matrix .
The first rows of (the right singular vectors) are returned in the array vt.
No rows of (no right singular vectors) are computed.
Constraint:
or .
3: – IntegerInput
On entry: , the number of rows of the matrix .
Constraint:
.
4: – IntegerInput
On entry: , the number of columns of the matrix .
Constraint:
.
5: – Real (Kind=nag_wp) arrayInput
Note: the second dimension of the array a
must be at least
.
On entry: the matrix .
6: – IntegerInput
On entry: the first dimension of the array a as declared in the (sub)program from which f10caf is called.
Constraint:
.
7: – IntegerInput
On entry: , number of columns in random projection, .
Constraint:
.
8: – Real (Kind=nag_wp)Input
On entry: the absolute tolerance, used in defining the threshold on estimating the rank of . If then is used unless in which case is used.
9: – Real (Kind=nag_wp)Input
On entry: the relative tolerance, used in defining the threshold on estimating the rank of . If then is used unless in which case is used.
10: – Integer arrayCommunication 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: – Real (Kind=nag_wp) arrayOutput
On exit: the first r elements of s contain the r largest singular values of in descending order. The remaining values are set to zero.
12: – Real (Kind=nag_wp) arrayOutput
Note: the second dimension of the array u
must be at least
if .
On exit: if , u contains the first r columns of (the left singular vectors, stored column-wise); the remaining elements of u are set to zero.
On entry: ifail must be set to , or to set behaviour on detection of an error; these values have no effect when no error is detected.
A value of causes the printing of an error message and program execution will be halted; otherwise program execution continues. A value of means that an error message is printed while a value of means that it is not.
If halting is not appropriate, the value or is recommended. If message printing is undesirable, then the value is recommended. Otherwise, the value is recommended. When the value or is used it is essential to test the value of ifail on exit.
On exit: unless the routine detects an error or a warning has been flagged (see Section 6).
6Error Indicators and Warnings
If on entry or , explanatory error messages are output on the current error message unit (as defined by x04aaf).
Errors or warnings detected by the routine:
On entry, . Constraint: or .
On entry, . Constraint: or .
On entry, . Constraint: .
On entry, . Constraint: .
On entry, . Constraint: .
On entry, . Constraint: .
On entry, state vector has been corrupted or not initialized.
On entry, and . Constraint: if , .
On entry, and . Constraint: if , .
On exit, , the rank of may be larger than r. Increase k to obtain a more accurate rank estimate. Smallest diagonal element of , from of , . Tolerance used to determine rank .
has effective rank of zero. First diagonal element of , from of , . Tolerance used to determine rank .
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.
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.
Dynamic memory allocation failed.
See Section 9 in the Introduction to the NAG Library FL Interface for further information.
7Accuracy
The error is approximately,
where,
The norm on the left-hand side of the first equation is the spectral norm, and is the th singular value of . More details on the error bound can be found in Sections 5 and 11 of Halko et al. (2011).
8Parallelism and Performance
Background information to multithreading can be found in the Multithreading documentation.
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.
9Further Comments
The total number of floating-point operations is . The first term corresponds to applying the random projection, i.e., computing . The second term corresponds to the decomposition of and the steps required to obtain the SVD of the original matrix .
Deterministic SVD solvers, such as f08kbf, require operations when and operations when .
The default values for rtol_abs and rtol_rel assume that you need an accurate approximation to . 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.
10Example
This example finds the singular values, the left and right singular vectors, and the rank of the matrix
using the randomised solver, f10caf, and a deterministic solver, f08kbf for comparison.