NAG CL Interface
f08khc (dgejsv)

1 Purpose

f08khc 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. f08khc implements the preconditioned Jacobi SVD of Drmač and Veselić (2008a) and Drmač and Veselić (2008b). This is the expert driver function that calls f08kjc after certain preconditioning. In most cases f08kbc or f08kdc, employing fast scaled rotations and de Rijk's pivoting strategy, is sufficient to obtain the SVD of a real matrix. These are much simpler to use and also handle the case m<n.

2 Specification

#include <nag.h>
void  f08khc (Nag_OrderType order, Nag_Preprocess joba, Nag_LeftVecsType jobu, Nag_RightVecsType jobv, Nag_ZeroCols jobr, Nag_TransType jobt, Nag_Perturb jobp, Integer m, Integer n, double a[], Integer pda, double sva[], double u[], Integer pdu, double v[], Integer pdv, double work[], Integer iwork[], NagError *fail)
The function may be called by the names: f08khc, nag_lapackeig_dgejsv or nag_dgejsv.

3 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, 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: order Nag_OrderType Input
On entry: the order argument specifies the two-dimensional storage scheme being used, i.e., row-major ordering or column-major ordering. C language defined storage is specified by order=Nag_RowMajor. See Section 3.1.3 in the Introduction to the NAG Library CL Interface for a more detailed explanation of the use of this argument.
Constraint: order=Nag_RowMajor or Nag_ColMajor.
2: joba Nag_Preprocess 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=Nag_ColpivRrank
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=Nag_ColpivRrankCond
Computation as with joba=Nag_ColpivRrank with an additional estimate of the condition number of B. It provides a realistic error bound.
joba=Nag_FullpivRrank
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=Nag_ColpivRrank option. If the structure of the input matrix is not known, and relative accuracy is desirable, then this option is advisable.
joba=Nag_FullpivRrankCond
Computation as with joba=Nag_FullpivRrank 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=Nag_ColpivSVrankAbs
Computation as with joba=Nag_ColpivRrank 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, is such that the relative residual norm (when comparing against A) is of the order Om×ε, where ε is machine precision. This gives the procedure licence to discard (set to zero) all singular values below n×ε×A.
joba=Nag_ColpivSVrankRel
Similar to joba=Nag_ColpivSVrankAbs. 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=Nag_ColpivSVrankAbs.
Constraint: joba=Nag_ColpivRrank, Nag_ColpivRrankCond, Nag_FullpivRrank, Nag_FullpivRrankCond, Nag_ColpivSVrankAbs or Nag_ColpivSVrankRel.
3: jobu Nag_LeftVecsType Input
On entry: specifies options for computing the left singular vectors U.
jobu=Nag_LeftSpan
The first n left singular vectors (columns of U) are computed and returned in the array u.
jobu=Nag_LeftVecs
All m left singular vectors are computed and returned in the array u.
jobu=Nag_NotLeftWork
No left singular vectors are computed, but the array u (with pdum and second dimension at least n) is available as workspace for computing right singular values. See the description of u.
jobu=Nag_NotLeftVecs
No left singular vectors are computed. u is not referenced when jobv=Nag_NotRightWork or Nag_NotRightVecs.
Constraint: jobu=Nag_LeftSpan, Nag_LeftVecs, Nag_NotLeftWork or Nag_NotLeftVecs.
4: jobv Nag_RightVecsType Input
On entry: specifies options for computing the right singular vectors V.
jobv=Nag_RightVecs
The n right singular vectors (columns of V) are computed and returned in the array v; Jacobi rotations are not explicitly accumulated.
jobv= Nag_RightVecsJRots
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=Nag_LeftSpan or Nag_LeftVecs, i.e., in computing the full SVD.
This is equivalent to multiplying the input matrix, on the right, by the matrix V.
jobv=Nag_NotRightWork
No right singular values are computed, but the array v (with pdvn and second dimension at least n) is available as workspace for computing left singular values. See the description of v.
jobv=Nag_NotRightVecs
No right singular vectors are computed. v is not referenced when jobu=Nag_NotLeftWork or Nag_NotLeftVecs or jobt=Nag_NoTrans or mn.
Constraints:
  • jobv=Nag_RightVecs, Nag_RightVecsJRots, Nag_NotRightWork or Nag_NotRightVecs;
  • if jobu=Nag_NotLeftWork or Nag_NotLeftVecs, jobv Nag_RightVecsJRots.
5: jobr Nag_ZeroCols 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 function to kill columns of A whose norm in cA is less than sfmin (for jobr=Nag_ZeroColsRestrict), or less than sfmin/ε (otherwise). sfmin is the safe range parameter, as returned by function X02AMC.
jobr=Nag_ZeroColsNormal
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 X02ALC, you are recommended to use function f08kjc.
jobr=Nag_ZeroColsRestrict
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 f08kjc.
Suggested value: jobr=Nag_ZeroColsRestrict.
Constraint: jobr=Nag_ZeroColsNormal or Nag_ZeroColsRestrict.
6: jobt Nag_TransType Input
On entry: 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 AT 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 ATA. See the descriptions of work[5] and work[6].
jobt=Nag_Trans
If n=m, perform an entropy test and, if the test indicates possibly faster convergence of the Jacobi process when using AT, then form the transpose AT. If A is replaced with AT, then the row pivoting is included automatically.
jobt=Nag_NoTrans
No entropy test and no transposition is performed.
The option jobt=Nag_Trans 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=Nag_Trans or Nag_NoTrans.
7: jobp Nag_Perturb Input
On entry: specifies whether the function 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=Nag_PerturbOn
Introduce perturbation if A is found to be very badly scaled (introducing denormalized numbers).
jobp=Nag_PerturbOff
Do not perturb.
Constraint: jobp=Nag_PerturbOn or Nag_PerturbOff.
8: m Integer Input
On entry: m, the number of rows of the matrix A.
Constraint: m0.
9: n Integer Input
On entry: n, the number of columns of the matrix A.
Constraint: mn0.
10: a[dim] double Input/Output
Note: the dimension, dim, of the array a must be at least
  • max1,pda×n when order=Nag_ColMajor;
  • max1,m×pda when order=Nag_RowMajor.
The i,jth element of the matrix A is stored in
  • a[j-1×pda+i-1] when order=Nag_ColMajor;
  • a[i-1×pda+j-1] when order=Nag_RowMajor.
On entry: the m by n matrix A.
On exit: the contents of a are overwritten.
11: pda Integer Input
On entry: the stride separating row or column elements (depending on the value of order) in the array a.
Constraints:
  • if order=Nag_ColMajor, pdamax1,m;
  • if order=Nag_RowMajor, pdamax1,n.
12: sva[n] double Output
On exit: the, possibly scaled, singular values of A.
The singular values of A are σi=α×sva[i-1], for i=1,2,,n, where α=work[0]/work[1]. 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=Nag_ZeroColsRestrict, then some of the singular values may be returned as exact zeros because they are below the numerical rank threshold or are denormalized numbers.
13: u[dim] double Output
Note: the dimension, dim, of the array u must be at least
  • max1,pdu×m when jobu=Nag_LeftVecs;
  • max1,pdu×n when jobu=Nag_LeftSpan or Nag_NotLeftWork and order=Nag_ColMajor;
  • max1,m×pdu when jobu=Nag_LeftSpan or Nag_NotLeftWork and order=Nag_RowMajor;
  • max1,m otherwise.
The i,jth element of the matrix U is stored in
  • u[j-1×pdu+i-1] when order=Nag_ColMajor;
  • u[i-1×pdu+j-1] when order=Nag_RowMajor.
On exit: if jobu=Nag_LeftSpan, u contains the m by n matrix of left singular vectors.
If jobu=Nag_LeftVecs, u contains the m by m matrix of left singular vectors, including an orthonormal basis of the orthogonal complement of Range(A).
u is not referenced when jobu=Nag_NotLeftWork or Nag_NotLeftVecs and one of the following is satisfied:
  • jobv=Nag_NotRightWork or Nag_NotRightVecs, or
  • n=1, or
  • A is the zero matrix.
14: pdu Integer Input
On entry: the stride separating row or column elements (depending on the value of order) in the array u.
Constraints:
  • if order=Nag_ColMajor,
    • if jobu=Nag_LeftVecs, pdumax1,m;
    • if jobu=Nag_LeftSpan or Nag_NotLeftWork, pdumax1,m;
    • otherwise pdu1;
  • if order=Nag_RowMajor,
    • if jobu=Nag_LeftVecs, pdumax1,m;
    • if jobu=Nag_LeftSpan or Nag_NotLeftWork, pdumax1,n;
    • otherwise pdu1.
15: v[dim] double Output
Note: the dimension, dim, of the array v must be at least
  • max1,pdv×n when jobv=Nag_RightVecs, Nag_RightVecsJRots or Nag_NotRightWork;
  • 1 otherwise.
The i,jth element of the matrix V is stored in
  • v[j-1×pdv+i-1] when order=Nag_ColMajor;
  • v[i-1×pdv+j-1] when order=Nag_RowMajor.
On exit: if jobv=Nag_RightVecs or Nag_RightVecsJRots, v contains the n by n matrix of right singular vectors.
v is not referenced when jobv=Nag_NotRightWork or Nag_NotRightVecs and one of the following is satisfied:
  • jobu=Nag_LeftSpan or Nag_LeftVecs and jobt=Nag_Trans, or
  • n=1, or
  • A is the zero matrix.
16: pdv Integer Input
On entry: the stride separating row or column elements (depending on the value of order) in the array v.
Constraints:
  • if jobv=Nag_RightVecs, Nag_RightVecsJRots or Nag_NotRightWork, pdvmax1,n;
  • otherwise pdv1.
17: work[7] double Output
On exit: contains information about the completed job.
work[0]
α=work[0]/work[1] is the scaling factor such that σi=α×sva[i-1], for i=1,2,,n are the computed singular values of A. (See the description of sva.)
work[1]
See the description of work[0].
work[2]
sconda, an estimate for the condition number of column equilibrated A (if joba=Nag_ColpivRrankCond or Nag_FullpivRrankCond). sconda is an estimate of RTR-11. It is computed using f07fgc. 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.
work[3]
An estimate of the scaled condition number of the triangular factor in the first QR factorization.
work[4]
An estimate of the scaled condition number of the triangular factor in the second QR factorization.
The following two parameters are computed if jobt=Nag_Trans.
work[5]
The entropy of ATA: this is the Shannon entropy of diagATA/traceATA taken as a point in the probability simplex.
work[6]
The entropy of AAT.
18: iwork[dim] Integer Output
On exit: contains information about the completed job.
iwork[0]
The numerical rank of A determined after the initial QR factorization with pivoting. See the descriptions of joba and jobr.
iwork[1]
The number of computed nonzero singular values.
iwork[2]
If nonzero, a warning message. If iwork[2]=1, then some of the column norms of A were denormalized (tiny) numbers. The requested high accuracy is not warranted by the data.
19: fail NagError * Input/Output
The NAG error argument (see Section 7 in the Introduction to the NAG Library CL Interface).

6 Error Indicators and Warnings

NE_ALLOC_FAIL
Dynamic memory allocation failed.
See Section 3.1.2 in the Introduction to the NAG Library CL Interface for further information.
NE_BAD_PARAM
On entry, argument value had an illegal value.
NE_CONSTRAINT
Constraint: jobv=Nag_RightVecs, Nag_RightVecsJRots, Nag_NotRightWork or Nag_NotRightVecs and
if jobu=Nag_NotLeftWork or Nag_NotLeftVecs, jobv Nag_RightVecsJRots.
NE_CONVERGENCE
f08khc did not converge in the allowed number of iterations (30). The computed values might be inaccurate.
NE_ENUM_INT_2
On entry, jobu=value, m=value and pdu=value.
Constraint: if jobu=Nag_LeftVecs, pdumax1,m;
if jobu=Nag_LeftSpan or Nag_NotLeftWork, pdumax1,m;
otherwise pdu1.
On entry, jobv=value, pdv=value and n=value.
Constraint: if jobv=Nag_RightVecs, Nag_RightVecsJRots or Nag_NotRightWork, pdvmax1,n;
otherwise pdv1.
NE_ENUM_INT_3
On entry, jobu=value, pdu=value, m=value and n=value.
Constraint: if jobu=Nag_LeftVecs, pdumax1,m;
if jobu=Nag_LeftSpan or Nag_NotLeftWork, pdumax1,n;
otherwise pdu1.
NE_INT
On entry, m=value.
Constraint: m0.
On entry, pda=value.
Constraint: pda>0.
On entry, pdu=value.
Constraint: pdu>0.
On entry, pdv=value.
Constraint: pdv>0.
NE_INT_2
On entry, m=value and n=value.
Constraint: mn0.
On entry, pda=value and m=value.
Constraint: pdamax1,m.
On entry, pda=value and n=value.
Constraint: pdamax1,n.
NE_INTERNAL_ERROR
An internal error has occurred in this function. Check the function call and any array sizes. If the call is correct then please contact NAG for assistance.
See Section 7.5 in the Introduction to the NAG Library CL Interface for further information.
NE_NO_LICENCE
Your licence key may have expired or may not have been installed correctly.
See Section 8 in the Introduction to the NAG Library CL Interface for further information.

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 orthogonal to working precision. See Section 4.9 of Anderson et al. (1999) for further details.

8 Parallelism and Performance

f08khc is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08khc 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 function. Please also consult the Users' Note for your implementation for any additional implementation-specific information.

9 Further Comments

f08khc implements a preconditioned Jacobi SVD algorithm. It uses f08aec, f08ahc and f08bfc 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 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.

10.1 Program Text

Program Text (f08khce.c)

10.2 Program Data

Program Data (f08khce.d)

10.3 Program Results

Program Results (f08khce.r)