F08KRF (ZGESDD) (PDF version)
F08 Chapter Contents
F08 Chapter Introduction
NAG Library Manual

NAG Library Routine Document

F08KRF (ZGESDD)

Note:  before using this routine, please read the Users' Note for your implementation to check the interpretation of bold italicised terms and other implementation-dependent details.

 Contents

    1  Purpose
    7  Accuracy

1  Purpose

F08KRF (ZGESDD) computes the singular value decomposition (SVD) of a complex m by n matrix A, 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
A = UΣVH ,  
where Σ is an m by n matrix which is zero except for its minm,n diagonal elements, U is an m by m unitary matrix, and V is an n by n unitary matrix. The diagonal elements of Σ are the singular values of A; they are real and non-negative, and are returned in descending order. The first minm,n columns of U and V are the left and right singular vectors of A.
Note that the routine returns VH, not V.

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:     JOBZ – CHARACTER(1)Input
On entry: specifies options for computing all or part of the matrix U.
JOBZ='A'
All m columns of U and all n rows of VH are returned in the arrays U and VT.
JOBZ='S'
The first minm,n columns of U and the first minm,n rows of VH are returned in the arrays U and VT.
JOBZ='O'
If MN, the first n columns of U are overwritten on the array A and all rows of VH are returned in the array VT. Otherwise, all columns of U are returned in the array U and the first m rows of VH are overwritten in the array VT.
JOBZ='N'
No columns of U or rows of VH are computed.
Constraint: JOBZ='A', 'S', 'O' or 'N'.
2:     M – INTEGERInput
On entry: m, the number of rows of the matrix A.
Constraint: M0.
3:     N – INTEGERInput
On entry: n, the number of columns of the matrix A.
Constraint: N0.
4:     ALDA* – COMPLEX (KIND=nag_wp) arrayInput/Output
Note: the second dimension of the array A must be at least max1,N.
On entry: the m by n matrix A.
On exit: if JOBZ='O', A is overwritten with the first n columns of U (the left singular vectors, stored column-wise) if MN; A is overwritten with the first m rows of VH (the right singular vectors, stored row-wise) otherwise.
If JOBZ'O', the contents of A are destroyed.
5:     LDA – INTEGERInput
On entry: the first dimension of the array A as declared in the (sub)program from which F08KRF (ZGESDD) is called.
Constraint: LDAmax1,M.
6:     SminM,N – REAL (KIND=nag_wp) arrayOutput
On exit: the singular values of A, sorted so that SiSi+1.
7:     ULDU* – COMPLEX (KIND=nag_wp) arrayOutput
Note: the second dimension of the array U must be at least max1,M if JOBZ='A' or JOBZ='O' and M<N, max1,minM,N if JOBZ='S', and at least 1 otherwise.
On exit:
If JOBZ='A' or JOBZ='O' and M<N, U contains the m by m unitary matrix U.
If JOBZ='S', U contains the first minm,n columns of U (the left singular vectors, stored column-wise).
If JOBZ='O' and MN, or JOBZ='N', U is not referenced.
8:     LDU – INTEGERInput
On entry: the first dimension of the array U as declared in the (sub)program from which F08KRF (ZGESDD) is called.
Constraints:
  • if JOBZ='S' or 'A' or JOBZ='O' and M<N, LDU max1,M ;
  • otherwise LDU1.
9:     VTLDVT* – COMPLEX (KIND=nag_wp) arrayOutput
Note: the second dimension of the array VT must be at least max1,N if JOBZ='A' or 'S' or JOBZ='O' and MN, and at least 1 otherwise.
On exit: if JOBZ='A' or JOBZ='O' and MN, VT contains the n by n unitary matrix VH.
If JOBZ='S', VT contains the first minm,n rows of VH (the right singular vectors, stored row-wise).
If JOBZ='O' and M<N, or JOBZ='N', VT is not referenced.
10:   LDVT – INTEGERInput
On entry: the first dimension of the array VT as declared in the (sub)program from which F08KRF (ZGESDD) is called.
Constraints:
  • if JOBZ='A' or JOBZ='O' and MN, LDVT max1,N ;
  • if JOBZ='S', LDVT max1,minM,N ;
  • otherwise LDVT1.
11:   WORKmax1,LWORK – COMPLEX (KIND=nag_wp) arrayWorkspace
On exit: if INFO=0, the real part of WORK1 contains the minimum value of LWORK required for optimal performance.
12:   LWORK – INTEGERInput
On entry: the dimension of the array WORK as declared in the (sub)program from which F08KRF (ZGESDD) is called.
If LWORK=-1, 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 nb×minM,N , where nb  is the optimal block size.
Constraints:
  • if JOBZ='N', LWORK2× minM,N+ max1,M,N ;
  • if JOBZ='O', LWORK2× minM,N× minM,N+2× minM,N+ max1,M,N ;
  • if JOBZ='S' or 'A', LWORK minM,N× minM,N+2× minM,N+ max1,M,N;
  • otherwise LWORK1.
13:   RWORK* – REAL (KIND=nag_wp) arrayWorkspace
Note: the dimension of the array RWORK must be at least max1,5×minM,N if JOBZ='N', and at least max1, minM,N× max 5× minM,N+7 , 2× maxM,N+ 2× minM,N+ 1  otherwise.
14:   IWORK8×minM,N – INTEGER arrayWorkspace
15:   INFO – INTEGEROutput
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
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 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

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.

9  Further Comments

The total number of floating-point operations is approximately proportional to mn2  when m>n and m2n  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 4 by 6 matrix
A = 0.96+0.81i -0.98-1.98i 0.62+0.46i -0.37-0.38i 0.83-0.51i 1.08+0.28i -0.03-0.96i -1.20-0.19i 1.01-0.02i 0.19+0.54i 0.20-0.01i 0.20+0.12i -0.91-2.06i -0.66-0.42i 0.63+0.17i -0.98+0.36i -0.17+0.46i -0.07-1.23i -0.05-0.41i -0.81-0.56i -1.11-0.60i 0.22+0.20i 1.47-1.59i 0.26-0.26i ,  
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 mn.

10.1  Program Text

Program Text (f08krfe.f90)

10.2  Program Data

Program Data (f08krfe.d)

10.3  Program Results

Program Results (f08krfe.r)


F08KRF (ZGESDD) (PDF version)
F08 Chapter Contents
F08 Chapter Introduction
NAG Library Manual

© The Numerical Algorithms Group Ltd, Oxford, UK. 2015