NAG FL Interface
f08kvf (zgejsv)

Settings help

FL Name Style:


FL Specification Language:


1 Purpose

f08kvf computes the singular value decomposition (SVD) of a real m×n matrix A, where mn, and optionally computes the left and/or right singular vectors. f08kvf implements the preconditioned Jacobi SVD of Drmač and Veselić (2008a) and Drmač and Veselić (2008b). This is an expert version of the Jacobi SVD routine f08kwf that employs preprocessing and preconditioning to improve the accuracy of results for extremely ill-conditioned matrices. In most cases f08kwf, employing fast scaled rotations and de Rijk's pivoting strategy, is sufficient and is simpler to use. Also consider using f08kpf or f08krf which are much simpler to use and also handle the case m<n.

2 Specification

Fortran Interface
Subroutine f08kvf ( joba, jobu, jobv, jobr, jobt, jobp, m, n, a, lda, sva, u, ldu, v, ldv, cwork, lwork, rwork, lrwork, iwork, info)
Integer, Intent (In) :: m, n, lda, ldu, ldv, lwork, lrwork
Integer, Intent (Out) :: iwork(m+3*n), info
Real (Kind=nag_wp), Intent (Out) :: sva(n), rwork(max(1,lrwork))
Complex (Kind=nag_wp), Intent (Inout) :: a(lda,*), u(ldu,*), v(ldv,*)
Complex (Kind=nag_wp), Intent (Out) :: cwork(lwork)
Character (1), Intent (In) :: joba, jobu, jobv, jobr, jobt, jobp
C Header Interface
#include <nag.h>
void  f08kvf_ (const char *joba, const char *jobu, const char *jobv, const char *jobr, const char *jobt, const char *jobp, const Integer *m, const Integer *n, Complex a[], const Integer *lda, double sva[], Complex u[], const Integer *ldu, Complex v[], const Integer *ldv, Complex cwork[], const Integer *lwork, double rwork[], const Integer *lrwork, Integer iwork[], Integer *info, const Charlen length_joba, const Charlen length_jobu, const Charlen length_jobv, const Charlen length_jobr, const Charlen length_jobt, const Charlen length_jobp)
The routine may be called by the names f08kvf, nagf_lapackeig_zgejsv or its LAPACK name zgejsv.

3 Description

The SVD is written as
A = UΣVH ,  
where Σ is an m×n matrix which is zero except for its n diagonal elements, U is an m×m unitary matrix, and V is an n×n unitary 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, respectively.

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
Drmač Z and Veselić K (2008a) New fast and accurate Jacobi SVD Algorithm I SIAM J. Matrix Anal. Appl. 29 4
Drmač Z and Veselić 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

5 Arguments

1: joba Character(1) Input
On entry: 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 can be written in the form 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 attempts to achieve 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ΣVH, is such that the relative residual norm (when comparing against A) is of the order O(m)×ε, 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 Character(1) Input
On entry: 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 when jobv='W' or 'N'.
Constraint: jobu='U', 'F', 'W' or 'N'.
3: jobv Character(1) Input
On entry: 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.
This is equivalent to multiplying the input matrix, on the right, by the matrix V.
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 when jobu='W' or 'N' or jobt='N' or mn.
Constraints:
  • jobv='V', 'J', 'W' or 'N';
  • if jobu='W' or 'N', jobv'J'.
4: jobr Character(1) Input
On entry: 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 routine 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 parameter, as returned by routine x02amf.
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 x02alf, you are recommended to use routine f08kwf.
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 f08kwf.
Suggested value: jobr='R'.
Constraint: jobr='N' or 'R'.
5: jobt Character(1) Input
On entry: specifies, in the case n=m, whether the routine is permitted to use the conjugate transpose of A for improved efficiency. If the matrix is square, then the procedure may use AH if it 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 AHA. See the descriptions of rwork(6) and rwork(7).
jobt='T'
If n=m, perform an entropy test and, if the test indicates possibly faster convergence of the Jacobi process when using AH, then form the conjugate transpose AH. If A is replaced with AH, 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 Character(1) Input
On entry: specifies whether the routine should be allowed to introduce structured perturbations to drown denormalized numbers. For details see Drmač and Veselić (2008a) and Drmač and Veselić (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: m Integer Input
On entry: m, the number of rows of the matrix A.
Constraint: m0.
8: n Integer Input
On entry: n, the number of columns of the matrix A.
Constraint: mn0.
9: a(lda,*) Complex (Kind=nag_wp) array Input/Output
Note: the second dimension of the array a must be at least max(1,n).
On entry: the m×n matrix A.
On exit: the contents of a are overwritten.
10: lda Integer Input
On entry: the first dimension of the array a as declared in the (sub)program from which f08kvf is called.
Constraint: ldamax(1,m).
11: sva(n) Real (Kind=nag_wp) array Output
On exit: the, possibly scaled, singular values of A.
The singular values of A are σi=α×sva(i), for i=1,2,,n, where α=rwork(1)/rwork(2). 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.
12: u(ldu,*) Complex (Kind=nag_wp) array Output
Note: the second dimension of the array u must be at least max(1,m) if jobu='F', max(1,n) if jobu='U' or 'W', and at least 1 otherwise.
On exit: if jobu='U', u contains the m×n matrix of left singular vectors.
If jobu='F', u contains the m×m matrix of left singular vectors, including an orthonormal basis of the orthogonal complement of Range(A).
u is not referenced when jobu='W' or 'N' and one of the following is satisfied:
  • jobv='W' or 'N', or
  • n=1, or
  • A is the zero matrix.
13: ldu Integer Input
On entry: the first dimension of the array u as declared in the (sub)program from which f08kvf is called.
Constraints:
  • if jobu='U', 'F' or 'W', ldumax(1,m);
  • otherwise ldu1.
14: v(ldv,*) Complex (Kind=nag_wp) array Output
Note: the second dimension of the array v must be at least max(1,n) if jobv='V', 'J' or 'W', and at least 1 otherwise.
On exit: if jobv='V' or 'J', v contains the n×n matrix of right singular vectors.
v is not referenced when jobv='W' or 'N' and one of the following is satisfied:
  • jobu='U' or 'F' and jobt='T', or
  • n=1, or
  • A is the zero matrix.
15: ldv Integer Input
On entry: the first dimension of the array v as declared in the (sub)program from which f08kvf is called.
Constraints:
  • if jobv='V', 'J' or 'W', ldvmax(1,n);
  • otherwise ldv1.
16: cwork(lwork) Complex (Kind=nag_wp) array Workspace
On exit: if info=0, cwork(1) contains the required minimal size of lwork.
17: lwork Integer Input
On entry: the dimension of the array cwork as declared in the (sub)program from which f08kvf is called.
If jobu='N' or 'W'
if jobv='N' or 'W'
if joba'E' or 'G'
the minimal requirement is lwork2n+1;
for optimal performance the requirement is lworkn+(n+1)×nb, where nb is the block size used by f08asf and f08btf. Assuming a value of nb=256 is wise, but choosing a smaller value (e.g., nb=128) should still lead to acceptable performance.
if joba='E' or 'G'
the minimal requirement is lworkn(n+3);
for optimal performance the requirement is lworkmax(n+(n+1)×nb,n(n+3)).
if jobv='V' or 'J'
the minimal requirement is lwork3n;
for optimal performance the requirement is lworkmax(nb,n)+n(nb+1).
If jobu='U' or 'F' and n1
if jobv='N' or 'W'
the minimal requirement is lwork3n;
for optimal performance the requirement is lworkmax(nb,n)+n(nb+1).
if jobv='V'
the minimal requirement is lworkn(2n+5);
for optimal performance the minimal requirement is sufficient.
if jobv='J'
the minimal requirement is lworkn(n+4);
for optimal performance the minimal requirement is sufficient.
If jobu='U' or 'F' and n=1
the minimal workspace requirement is lworkm+1;
for optimal performance, lworkm×nb+1.
18: rwork(max(1,lrwork)) Real (Kind=nag_wp) array Workspace
On exit: contains information about the completed job.
rwork(1)
α=rwork(1)/rwork(2) is the scaling factor such that σi=α×sva(i), for i=1,2,,n, are the computed singular values of A. See the description of sva.
rwork(2)
See the description of rwork(1).
rwork(3)
sconda, an estimate for the condition number of column equilibrated A (if joba='E' or 'G'). sconda is an estimate of ((RHR)-11). It is computed using f07fgf. 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 the 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.
rwork(4)
An estimate of the scaled condition number of the triangular factor in the first QR factorization.
rwork(5)
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'.
rwork(6)
The entropy of AHA: this is the Shannon entropy of diagAHA/traceAHA taken as a point in the probability simplex.
rwork(7)
The entropy of AAH.
19: lrwork Integer Input
On entry: the first dimension of the array rwork as declared in the (sub)program from which f08kvf is called.
Constraint: lrworkmax(7,n+2×m).
20: iwork(m+3×n) Integer array Workspace
On exit: contains information about the completed job.
iwork(1)
The numerical rank of A determined after the initial QR factorization with pivoting. See the descriptions of joba and jobr.
iwork(2)
The number of computed nonzero singular values.
iwork(3)
If nonzero, a warning message. If iwork(3)=1, then some of the column norms of A were denormalized (tiny) numbers. The requested high accuracy is not warranted by the data.
21: info Integer Output
On exit: info=0 unless the routine detects an error (see Section 6).

6 Error Indicators and Warnings

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

7 Accuracy

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

8 Parallelism and Performance

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

f08kvf implements a preconditioned Jacobi SVD algorithm. It uses f08asf, f08avf and f08btf 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 a 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.

10 Example

This example finds the singular values and left and right singular vectors of the 6×4 matrix
A = ( 0.96-0.81i -0.03+0.96i -0.91+2.06i -0.05+0.41i -0.98+1.98i -1.20+0.19i -0.66+0.42i -0.81+0.56i 0.62-0.46i 1.01+0.02i 0.63-0.17i -1.11+0.60i -0.37+0.38i 0.19-0.54i -0.98-0.36i 0.22-0.20i 0.83+0.51i 0.20+0.01i -0.17-0.46i 1.47+1.59i 1.08-0.28i 0.20-0.12i -0.07+1.23i 0.26+0.26i ) ,  
together with the condition number of A and approximate error bounds for the computed singular values and vectors.

10.1 Program Text

Program Text (f08kvfe.f90)

10.2 Program Data

Program Data (f08kvfe.d)

10.3 Program Results

Program Results (f08kvfe.r)