?gesdd

Computes the singular value decomposition of a general rectangular matrix using a divide and conquer method.

Syntax

FORTRAN 77:

call sgesdd(jobz, m, n, a, lda, s, u, ldu, vt, ldvt, work, lwork, iwork, info)

call dgesdd(jobz, m, n, a, lda, s, u, ldu, vt, ldvt, work, lwork, iwork, info)

call cgesdd(jobz, m, n, a, lda, s, u, ldu, vt, ldvt, work, lwork, rwork, iwork, info)

call zgesdd(jobz, m, n, a, lda, s, u, ldu, vt, ldvt, work, lwork, rwork, iwork, info)

Fortran 95:

call gesdd(a, s [,u] [,vt] [,jobz] [,info])

C:

lapack_int LAPACKE_sgesdd( int matrix_order, char jobz, lapack_int m, lapack_int n, float* a, lapack_int lda, float* s, float* u, lapack_int ldu, float* vt, lapack_int ldvt );

lapack_int LAPACKE_dgesdd( int matrix_order, char jobz, lapack_int m, lapack_int n, double* a, lapack_int lda, double* s, double* u, lapack_int ldu, double* vt, lapack_int ldvt );

lapack_int LAPACKE_cgesdd( int matrix_order, char jobz, lapack_int m, lapack_int n, lapack_complex_float* a, lapack_int lda, float* s, lapack_complex_float* u, lapack_int ldu, lapack_complex_float* vt, lapack_int ldvt );

lapack_int LAPACKE_zgesdd( int matrix_order, char jobz, lapack_int m, lapack_int n, lapack_complex_double* a, lapack_int lda, double* s, lapack_complex_double* u, lapack_int ldu, lapack_complex_double* vt, lapack_int ldvt );

Include Files

The FORTRAN 77 interfaces are specified in the mkl_lapack.fi and mkl_lapack.h include files, the Fortran 95 interfaces are specified in the lapack.f90 include file, and the C interfaces are specified in the mkl_lapacke.h include file.

Description

The routine computes the singular value decomposition (SVD) of a real/complex m-by-n matrix A, optionally computing the left and/or right singular vectors.

If singular vectors are desired, it uses a divide-and-conquer algorithm. The SVD is written

A = U*Σ*V' for real routines,

A = U*Σ*conjg(V') for complex routines,

where Σ is an m-by-n matrix which is zero except for its min(m,n) diagonal elements, U is an m-by-m orthogonal/unitary matrix, and V is an n-by-n orthogonal/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 min(m, n) columns of U and V are the left and right singular vectors of A.

Note that the routine returns vt = V' (for real flavors) or vt = conjg(V') (for complex flavors), not V.

Input Parameters

The data types are given for the Fortran interface. A <datatype> placeholder, if present, is used for the C interface data types in the C interface section above. See C Interface Conventions for the C interface principal conventions and type definitions.

jobz

CHARACTER*1. Must be 'A', 'S', 'O', or 'N'.

Specifies options for computing all or part of the matrix U.

If jobz = 'A', all m columns of U and all n rows of V'/conjg(V') are returned in the arrays u and vt;

if jobz = 'S', the first min(m, n) columns of U and the first min(m, n) rows of V'/conjg(V') are returned in the arrays u and vt;

if jobz = 'O', then

if m   n, the first n columns of U are overwritten in the array a and all rows of V'/conjg(V') are returned in the array vt;

if m < n, all columns of U are returned in the array u and the first m rows of V'/conjg(V') are overwritten in the array a;

if jobz = 'N', no columns of U or rows of V'/conjg(V') are computed.

m

INTEGER. The number of rows of the matrix A (m 0).

n

INTEGER. The number of columns in A (n  0).

a, work

REAL for sgesdd

DOUBLE PRECISION for dgesdd

COMPLEX for cgesdd

DOUBLE COMPLEX for zgesdd.

Arrays: a(lda,*) is an array containing the m-by-n matrix A.

The second dimension of a must be at least max(1, n).

work is a workspace array, its dimension max(1, lwork).

lda

INTEGER. The leading dimension of the array a. Must be at least max(1, m).

ldu, ldvt

INTEGER. The leading dimensions of the output arrays u and vt, respectively.

Constraints:

ldu 1; ldvt 1.

If jobz = 'S' or 'A', or jobz = 'O' and m < n,

then ldu m;

If jobz = 'A' or jobz = 'O' and m  n,

then ldvt n;

If jobz = 'S', ldvt min(m, n).

lwork

INTEGER.

The dimension of the array work; lwork 1.

If lwork = -1, then a workspace query is assumed; the routine only calculates the optimal size of the work array, returns this value as the work(1), and no error message related to lwork is issued by xerbla.

See Application Notes for the suggested value of lwork.

rwork

REAL for cgesdd

DOUBLE PRECISION for zgesdd

Workspace array, DIMENSION at least max(1, 5*min(m,n)) if jobz = 'N'.

Otherwise, the dimension of rwork must be at least max(1,min(m,n)*max(5*min(m,n)+7,2*max(m,n)+2*min(m,n)+1)). This array is used in complex flavors only.

iwork

INTEGER. Workspace array, DIMENSION at least max(1, 8 *min(m, n)).

Output Parameters

a

On exit:

If jobz = 'O', then if m n, a is overwritten with the first n columns of U (the left singular vectors, stored columnwise). If m < n, a is overwritten with the first m rows of VT (the right singular vectors, stored rowwise);

If jobz'O', the contents of a are destroyed.

s

REAL for single precision flavors DOUBLE PRECISION for double precision flavors.

Array, DIMENSION at least max(1, min(m,n)). Contains the singular values of A sorted so that s(i)  s(i+1).

u, vt

REAL for sgesdd

DOUBLE PRECISION for dgesdd

COMPLEX for cgesdd

DOUBLE COMPLEX for zgesdd.

Arrays:

u(ldu,*); the second dimension of u must be at least max(1, m) if jobz = 'A' or jobz = 'O' and m < n.

If jobz = 'S', the second dimension of u must be at least max(1, min(m, n)).

If jobz = 'A'or jobz = 'O' and m < n, u contains the m-by-m orthogonal/unitary matrix U.

If jobz = 'S', u contains the first min(m, n) columns of U (the left singular vectors, stored columnwise).

If jobz = 'O' and mn, or jobz = 'N', u is not referenced.

vt(ldvt,*); the second dimension of vt must be at least max(1, n).

If jobz = 'A'or jobz = 'O' and mn, vt contains the n-by-n orthogonal/unitary matrix VT.

If jobz = 'S', vt contains the first min(m, n) rows of VT (the right singular vectors, stored rowwise).

If jobz = 'O' and m < n, or jobz = 'N', vt is not referenced.

work(1)

On exit, if info = 0, then work(1) returns the required minimal size of lwork.

info

INTEGER.

If info = 0, the execution is successful.

If info = -i, the i-th parameter had an illegal value.

If info = i, then ?bdsdc did not converge, updating process failed.

Fortran 95 Interface Notes

Routines in Fortran 95 interface have fewer arguments in the calling sequence than their FORTRAN 77 counterparts. For general conventions applied to skip redundant or restorable arguments, see Fortran 95 Interface Conventions.

Specific details for the routine gesdd interface are the following:

a

Holds the matrix A of size (m, n).

s

Holds the vector of length min(m, n).

u

Holds the matrix U of size

  • (m,m) if jobz='A' or jobz='O' and m < n

  • (m,min(m, n)) if jobz='S'

u is not referenced if jobz is not supplied or if jobz='N' or jobz='O' and m n.

vt

Holds the matrix VT of size

  • (n,n) if jobz='A' or jobz='O' and m n

  • (min(m, n), n) if jobz='S'

vt is not referenced if jobz is not supplied or if jobz='N' or jobz='O' and m < n.

job

Must be 'N', 'A', 'S', or 'O'. The default value is 'N'.

Application Notes

For real flavors:

If jobz = 'N', lwork 3*min(m, n) + max (max(m,n), 6*min(m, n));

If jobz = 'O', lwork 3*(min(m, n))2 + max (max(m, n), 5*(min(m, n))2 + 4*min(m, n));

If jobz = 'S' or 'A', lwork 3*(min(m, n))2 + max (max(m, n), 4*(min(m, n))2 + 4*min(m, n))

For complex flavors:

If jobz = 'N', lwork 2*min(m, n) + max(m, n);

If jobz = 'O', lwork 2*(min(m, n))2 + max(m, n) + 2*min(m, n);

If jobz = 'S' or 'A', lwork (min(m, n))2 + max(m, n) + 2*min(m, n);

For good performance, lwork should generally be larger.

If you are in doubt how much workspace to supply, use a generous value of lwork for the first run or set lwork = -1.

If you choose the first option and set any of admissible lwork sizes, which is no less than the minimal value described, the routine completes the task, though probably not so fast as with a recommended workspace, and provides the recommended workspace in the first element of the corresponding array work on exit. Use this value (work(1)) for subsequent runs.

If you set lwork = -1, the routine returns immediately and provides the recommended workspace in the first element of the corresponding array (work). This operation is called a workspace query.

Note that if you set lwork to less than the minimal required value and not -1, the routine returns immediately with an error exit and does not provide any information on the recommended workspace.


Submit feedback on this help topic

Copyright © 1994 - 2011, Intel Corporation. All rights reserved.