NAG FL Interface
f08khf (dgejsv)
1
Purpose
f08khf computes the singular value decomposition (SVD) of a real
by
matrix
, where
, and optionally computes the left and/or right singular vectors.
f08khf implements the preconditioned Jacobi SVD of
Drmač and Veselić (2008a) and
Drmač and Veselić (2008b). This is the expert driver routine that calls
f08kjf after certain preconditioning. In most cases
f08kbf or
f08kdf, 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
.
2
Specification
Fortran Interface
Subroutine f08khf ( |
joba, jobu, jobv, jobr, jobt, jobp, m, n, a, lda, sva, u, ldu, v, ldv, work, lwork, iwork, info) |
Integer, Intent (In) |
:: |
m, n, lda, ldu, ldv, lwork |
Integer, Intent (Out) |
:: |
iwork(m+3*n), info |
Real (Kind=nag_wp), Intent (Inout) |
:: |
a(lda,*), u(ldu,*), v(ldv,*) |
Real (Kind=nag_wp), Intent (Out) |
:: |
sva(n), work(lwork) |
Character (1), Intent (In) |
:: |
joba, jobu, jobv, jobr, jobt, jobp |
|
C Header Interface
#include <nag.h>
void |
f08khf_ (const char *joba, const char *jobu, const char *jobv, const char *jobr, const char *jobt, const char *jobp, const Integer *m, const Integer *n, double a[], const Integer *lda, double sva[], double u[], const Integer *ldu, double v[], const Integer *ldv, double work[], const Integer *lwork, 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) |
|
C++ Header Interface
#include <nag.h> extern "C" {
void |
f08khf_ (const char *joba, const char *jobu, const char *jobv, const char *jobr, const char *jobt, const char *jobp, const Integer &m, const Integer &n, double a[], const Integer &lda, double sva[], double u[], const Integer &ldu, double v[], const Integer &ldv, double work[], const Integer &lwork, 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 f08khf, nagf_lapackeig_dgejsv or its LAPACK name dgejsv.
3
Description
The SVD is written as
where
is an
by
matrix which is zero except for its
diagonal elements,
is an
by
orthogonal matrix, and
is an
by
orthogonal matrix. The diagonal elements of
are the singular values of
in descending order of magnitude. The columns of
and
are the left and the right singular vectors of
, 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:
– Character(1)
Input
-
On entry: specifies the form of pivoting for the
factorization stage; whether an estimate of the condition number of the scaled matrix is required; and the form of rank reduction that is performed.
- The initial 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 . This option works well (high relative accuracy) if can be written in the form , with well-conditioned and arbitrary diagonal matrix . The accuracy cannot be spoiled by column scaling. The accuracy of the computed output depends on the condition of , and the procedure attempts to achieve the best theoretical accuracy.
- Computation as with with an additional estimate of the condition number of . It provides a realistic error bound.
- The initial 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 . If with ill-conditioned diagonal scalings , , and well-conditioned matrix , this option gives higher accuracy than the option. If the structure of the input matrix is not known, and relative accuracy is desirable, then this option is advisable.
- Computation as with with an additional estimate of the condition number of , where (i.e., ). If has heavily weighted rows, then using this condition number gives too pessimistic an error bound.
- Computation as with 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, , is such that the relative residual norm (when comparing against ) is of the order , where is machine precision. This gives the procedure licence to discard (set to zero) all singular values below .
- Similar to . The rank revealing property of the initial factorization is used to reveal (using the upper triangular factor) a gap, , in which case the numerical rank is declared to be . The SVD is computed with absolute error bounds, but more accurately than with .
Constraint:
, , , , or .
-
2:
– Character(1)
Input
-
On entry: specifies options for computing the left singular vectors
.
- The first left singular vectors (columns of ) are computed and returned in the array u.
- All left singular vectors are computed and returned in the array u.
- No left singular vectors are computed, but the array u (with and second dimension at least n) is available as workspace for computing right singular values. See the description of u.
- No left singular vectors are computed. is not referenced when or .
Constraint:
, , or .
-
3:
– Character(1)
Input
-
On entry: specifies options for computing the right singular vectors
.
- The right singular vectors (columns of ) are computed and returned in the array v; Jacobi rotations are not explicitly accumulated.
- The right singular vectors (columns of ) are computed and returned in the array v, but they are computed as the product of Jacobi rotations. This option is allowed only if or , i.e., in computing the full SVD.
This is equivalent to multiplying the input matrix, on the right, by the matrix .
- No right singular values are computed, but the array v (with and second dimension at least n) is available as workspace for computing left singular values. See the description of v.
- No right singular vectors are computed. is not referenced when or or or .
Constraints:
- , , or ;
- if or , .
-
4:
– Character(1)
Input
-
On entry: specifies the conditions under which columns of
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
is scaled so that the largest column (in the Euclidean norm) of
is equal to the square root of the overflow threshold, then
jobr allows the routine to kill columns of
whose norm in
is less than
(for
), or less than
(otherwise).
is the safe range parameter, as returned by routine
x02amf.
- Only set to zero those columns of for which the norm of corresponding column of , that is, those columns that are effectively zero (to machine precision) anyway. If the condition number of is greater than the overflow threshold , where is the value returned by x02alf, you are recommended to use routine f08kjf.
- Set to zero those columns of for which the norm of the corresponding column of . This approximately represents a restricted range for of .
For computing the singular values in the full range from the safe minimum up to the overflow threshold use
f08kjf.
Suggested value:
.
Constraint:
or .
-
5:
– Character(1)
Input
-
On entry: specifies, in the case
, whether the routine is permitted to use the transpose of
for improved efficiency. If the matrix is square, then the procedure may use
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
. See the descriptions of
and
.
- If , perform an entropy test and, if the test indicates possibly faster convergence of the Jacobi process when using , then form the transpose . If is replaced with , then the row pivoting is included automatically.
- No entropy test and no transposition is performed.
The option
can be used to compute only the singular values, or the full SVD (
,
and
). In the case where only one set of singular vectors (
or
) is required, the caller must still provide both
u and
v, as one of the matrices is used as workspace if the matrix
is transposed. See the descriptions of
u and
v.
Constraint:
or .
-
6:
– 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.
- Introduce perturbation if is found to be very badly scaled (introducing denormalized numbers).
- Do not perturb.
Constraint:
or .
-
7:
– Integer
Input
-
On entry: , the number of rows of the matrix .
Constraint:
.
-
8:
– Integer
Input
-
On entry: , the number of columns of the matrix .
Constraint:
.
-
9:
– Real (Kind=nag_wp) array
Input/Output
-
Note: the second dimension of the array
a
must be at least
.
On entry: the by matrix .
On exit: the contents of
a are overwritten.
-
10:
– Integer
Input
-
On entry: the first dimension of the array
a as declared in the (sub)program from which
f08khf is called.
Constraint:
.
-
11:
– Real (Kind=nag_wp) array
Output
-
On exit: the, possibly scaled, singular values of
.
The singular values of are
, for , where . Normally and no scaling is required to obtain the singular values. However, if the largest singular value of overflows or if small singular values have been saved from underflow by scaling the input matrix , then .
If , 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:
– Real (Kind=nag_wp) array
Output
-
Note: the second dimension of the array
u
must be at least
if
,
if
or
, and at least
otherwise.
On exit: if
,
u contains the
by
matrix of left singular vectors.
If
,
u contains the
by
matrix of left singular vectors, including an orthonormal basis of the orthogonal complement of
Range().
u is not referenced when
or
and one of the following is satisfied:
- or , or
- , or
- is the zero matrix.
-
13:
– Integer
Input
-
On entry: the first dimension of the array
u as declared in the (sub)program from which
f08khf is called.
Constraints:
- if , or , ;
- otherwise .
-
14:
– Real (Kind=nag_wp) array
Output
-
Note: the second dimension of the array
v
must be at least
if
,
or
, and at least
otherwise.
On exit: if
or
,
v contains the
by
matrix of right singular vectors.
v is not referenced when
or
and one of the following is satisfied:
- or and , or
- , or
- is the zero matrix.
-
15:
– Integer
Input
-
On entry: the first dimension of the array
v as declared in the (sub)program from which
f08khf is called.
Constraints:
- if , or , ;
- otherwise .
-
16:
– Real (Kind=nag_wp) array
Workspace
-
On exit: contains information about the completed job.
- is the scaling factor such that
, for are the computed singular values of . (See the description of .)
- See the description of .
- sconda, an estimate for the condition number of column equilibrated (if or ). sconda is an estimate of . It is computed using f07fgf. It satisfies where is the triangular factor from the factorization of . However, if is truncated and the numerical rank is determined to be strictly smaller than , sconda is returned as , 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.
- An estimate of the scaled condition number of the triangular factor in the first factorization.
- An estimate of the scaled condition number of the triangular factor in the second factorization.
The following two parameters are computed if .
- The entropy of : this is the Shannon entropy of taken as a point in the probability simplex.
- The entropy of .
-
17:
– Integer
Input
-
On entry: the dimension of the array
work as declared in the (sub)program from which
f08khf is called.
- If or and or
-
- if or
-
the minimal requirement is ;
for optimal performance the requirement is
, where
is the block size used by
f08aef and
f08bff. Assuming a value of
is wise, but choosing a smaller value (e.g.,
) should still lead to acceptable performance.
- if or
-
the minimal requirement is ;
for optimal performance the requirement is .
- If or and or
-
the minimal requirement is ;
for optimal performance the requirement is , where is described above.
- If or and or
-
the minimal requirement is ;
for optimal performance the requirement is , where is described above.
- If or and
-
.
- If or and
-
the minimal requirement is ;
for better performance , where is described above.
-
18:
– Integer array
Workspace
-
On exit: contains information about the completed job.
- The numerical rank of determined after the initial factorization with pivoting. See the descriptions of joba and jobr.
- The number of computed nonzero singular values.
- If nonzero, a warning message: If , then some of the column norms of were denormalized (tiny) numbers. The requested high accuracy is not warranted by the data.
-
19:
– Integer
Output
On exit:
unless the routine detects an error (see
Section 6).
6
Error Indicators and Warnings
If , argument had an illegal value. An explanatory message is output, and execution of the program is terminated.
-
f08khf did not converge in the allowed number of iterations (). The computed values might be inaccurate.
7
Accuracy
The computed SVD is nearly the exact SVD for a nearby matrix
, where
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
f08khf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f08khf 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.
f08khf implements a preconditioned Jacobi SVD algorithm. It uses
f08aef,
f08ahf and
f08bff 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
with the structure
, where
,
are arbitrarily ill-conditioned diagonal matrices and
is a well-conditioned matrix. In that case, complete pivoting in the first
factorizations provides accuracy dependent on the condition number of
, and independent of
,
. Such higher accuracy is not completely understood theoretically, but it works well in practice. Further, if
can be written as
, with well-conditioned
and some diagonal
, then the high accuracy is guaranteed, both theoretically and in software, independent of
.
10
Example
This example finds the singular values and left and right singular vectors of the
by
matrix
together with the condition number of
and approximate error bounds for the computed singular values and vectors.
10.1
Program Text
10.2
Program Data
10.3
Program Results