dgesdd | R Documentation |
DGESDD computes the singular value decomposition (SVD) of a real M-by-N matrix A, optionally computing the left and right singular vectors. If singular vectors are desired, it uses a divide-and-conquer algorithm.
The SVD is written
A = U * SIGMA * transpose(V)
where SIGMA is an M-by-N matrix which is zero except for its min(m,n) diagonal elements, U is an M-by-M orthogonal matrix, and V is an N-by-N orthogonal matrix. The diagonal elements of SIGMA 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**T, not V.
dgesdd(
JOBZ = "A",
M = NULL,
N = NULL,
A,
LDA = NULL,
S,
U,
LDU = NULL,
VT,
LDVT = NULL,
WORK = NULL,
LWORK = NULL
)
JOBZ |
a character. Specifies options for computing all or part of the matrix U:
|
M |
an integer. The number of rows of the input matrix A. M >= 0. |
N |
an integer. The number of columns of the input matrix A. N >= 0. |
A |
the M-by-N matrix A. |
LDA |
an integer. The leading dimension of the matrix A. LDA >= max(1,M). |
S |
a matrix of dimension (min(M,N)). The singular values of A, sorted so that S(i) >= S(i+1). |
U |
U is a matrx of dimension (LDU,UCOL)
|
LDU |
an integer. The leading dimension of the matrix U. LDU >= 1; if JOBZ = 'S' or 'A' or JOBZ = 'O' and M < N, LDU >= M. |
VT |
VT is matrix of dimension (LDVT,N)
|
LDVT |
an integer. The leading dimension of the matrix VT. LDVT >= 1; if JOBZ = 'A' or JOBZ = 'O' and M >= N, LDVT >= N; if JOBZ = 'S', LDVT >= min(M,N). |
WORK |
a matrix of dimension (MAX(1,LWORK)) |
LWORK |
an integer. The dimension of the array WORK. LWORK >= 1. If LWORK = -1, a workspace query is assumed. The optimal size for the WORK array is calculated and stored in WORK(1), and no other work except argument checking is performed. Let mx = max(M,N) and mn = min(M,N).
These are not tight minimums in all cases; see comments inside code. For good performance, LWORK should generally be larger; a query is recommended. |
IWORK an integer matrix dimension of (8*min(M,N)) A is updated.
JOBZ = 'O', 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 V**T (the right singular vectors, stored rowwise) otherwise.
JOBZ .ne. 'O', the contents of A are destroyed.
INFO an integer
successful exit.
if INFO = -i, the i-th argument had an illegal value.
DBDSDC did not converge, updating process failed.
## Not run:
set.seed(4669)
A = matrix(rnorm(12),4,3)
S = matrix(0,nrow=3,ncol=1)
U = matrix(0,nrow=4,ncol=4)
VT = matrix(0,ncol=3,nrow=3)
dgesdd(A=A,S=S,U=U,VT=VT)
S
U
VT
rm(A,S,U,VT)
A = as.big.matrix(matrix(rnorm(12),4,3))
S = as.big.matrix(matrix(0,nrow=3,ncol=1))
U = as.big.matrix(matrix(0,nrow=4,ncol=4))
VT = as.big.matrix(matrix(0,ncol=3,nrow=3))
dgesdd(A=A,S=S,U=U,VT=VT)
S[,]
U[,]
VT[,]
rm(A,S,U,VT)
gc()
## End(Not run)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.