NAG Library Routine Document
F08KRF (ZGESDD)
1 Purpose
F08KRF (ZGESDD) computes the singular value decomposition (SVD) of a complex by matrix , optionally computing the left and/or right singular vectors, by using a divide-and-conquer method.
2 Specification
SUBROUTINE F08KRF ( |
JOBZ, M, N, A, LDA, S, U, LDU, VT, LDVT, WORK, LWORK, RWORK, IWORK, INFO) |
INTEGER |
M, N, LDA, LDU, LDVT, LWORK, IWORK(8*min(M,N)), INFO |
REAL (KIND=nag_wp) |
S(min(M,N)), RWORK(*) |
COMPLEX (KIND=nag_wp) |
A(LDA,*), U(LDU,*), VT(LDVT,*), WORK(max(1,LWORK)) |
CHARACTER(1) |
JOBZ |
|
The routine may be called by its
LAPACK
name zgesdd.
3 Description
The SVD is written as
where
is an
by
matrix which is zero except for its
diagonal elements,
is an
by
unitary matrix, and
is an
by
unitary 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 .
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
http://www.netlib.org/lapack/lug
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore
5 Parameters
- 1: – CHARACTER(1)Input
-
On entry: specifies options for computing all or part of the matrix
.
- All columns of and all rows of are returned in the arrays U and VT.
- The first columns of and the first rows of are returned in the arrays U and VT.
- If , the first columns of are overwritten on the array A and all rows of are returned in the array VT. Otherwise, all columns of are returned in the array U and the first rows of are overwritten in the array VT.
- No columns of or rows of are computed.
Constraint:
, , or .
- 2: – INTEGERInput
-
On entry: , the number of rows of the matrix .
Constraint:
.
- 3: – INTEGERInput
-
On entry: , the number of columns of the matrix .
Constraint:
.
- 4: – COMPLEX (KIND=nag_wp) arrayInput/Output
-
Note: the second dimension of the array
A
must be at least
.
On entry: the by matrix .
On exit: if
,
A is overwritten with the first
columns of
(the left singular vectors, stored column-wise) if
;
A is overwritten with the first
rows of
(the right singular vectors, stored row-wise) otherwise.
If
, the contents of
A are destroyed.
- 5: – INTEGERInput
-
On entry: the first dimension of the array
A as declared in the (sub)program from which F08KRF (ZGESDD) is called.
Constraint:
.
- 6: – REAL (KIND=nag_wp) arrayOutput
-
On exit: the singular values of , sorted so that .
- 7: – COMPLEX (KIND=nag_wp) arrayOutput
-
Note: the second dimension of the array
U
must be at least
if
or
and
,
if
, and at least
otherwise.
On exit:
If
or
and
,
U contains the
by
unitary matrix
.
If
,
U contains the first
columns of
(the left singular vectors, stored column-wise).
If
and
, or
,
U is not referenced.
- 8: – INTEGERInput
-
On entry: the first dimension of the array
U as declared in the (sub)program from which F08KRF (ZGESDD) is called.
Constraints:
- if or or and , ;
- otherwise .
- 9: – COMPLEX (KIND=nag_wp) arrayOutput
-
Note: the second dimension of the array
VT
must be at least
if
or
or
and
, and at least
otherwise.
On exit: if
or
and
,
VT contains the
by
unitary matrix
.
If
,
VT contains the first
rows of
(the right singular vectors, stored row-wise).
If
and
, or
,
VT is not referenced.
- 10: – INTEGERInput
-
On entry: the first dimension of the array
VT as declared in the (sub)program from which F08KRF (ZGESDD) is called.
Constraints:
- if or and , ;
- if , ;
- otherwise .
- 11: – COMPLEX (KIND=nag_wp) arrayWorkspace
-
On exit: if
, the real part of
contains the minimum value of
LWORK required for optimal performance.
- 12: – INTEGERInput
-
On entry: the dimension of the array
WORK as declared in the (sub)program from which F08KRF (ZGESDD) is called.
If
, a workspace query is assumed; the routine only calculates the optimal size of the
WORK array, returns this value as the first entry of the
WORK array, and no error message related to
LWORK is issued.
Suggested value:
for optimal performance,
LWORK should generally be larger. Consider increasing
LWORK by at least
, where
is the optimal
block size.
Constraints:
- if , ;
- if , ;
- if or , ;
- otherwise .
- 13: – REAL (KIND=nag_wp) arrayWorkspace
-
Note: the dimension of the array
RWORK
must be at least
if
, and at least
otherwise.
- 14: – INTEGER arrayWorkspace
-
- 15: – INTEGEROutput
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.
-
F08KRF (ZGESDD) did not converge, the updating process failed.
7 Accuracy
The computed singular value decomposition is nearly the exact singular value decomposition 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
F08KRF (ZGESDD) is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
F08KRF (ZGESDD) 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.
The total number of floating-point operations is approximately proportional to when and otherwise.
The singular values are returned in descending order.
The real analogue of this routine is
F08KDF (DGESDD).
10 Example
This example finds the singular values and left and right singular vectors of the
by
matrix
together with approximate error bounds for the computed singular values and vectors.
The example program for
F08KPF (ZGESVD) illustrates finding a singular value decomposition for the case
.
10.1 Program Text
Program Text (f08krfe.f90)
10.2 Program Data
Program Data (f08krfe.d)
10.3 Program Results
Program Results (f08krfe.r)