Skip to main content

cusolverDnSgesvdaStridedBatched

Function cusolverDnSgesvdaStridedBatched 

Source
pub unsafe extern "C" fn cusolverDnSgesvdaStridedBatched(
    handle: cusolverDnHandle_t,
    jobz: cusolverEigMode_t,
    rank: c_int,
    m: c_int,
    n: c_int,
    d_A: *const f32,
    lda: c_int,
    strideA: c_longlong,
    d_S: *mut f32,
    strideS: c_longlong,
    d_U: *mut f32,
    ldu: c_int,
    strideU: c_longlong,
    d_V: *mut f32,
    ldv: c_int,
    strideV: c_longlong,
    d_work: *mut f32,
    lwork: c_int,
    d_info: *mut c_int,
    h_R_nrmF: *mut f64,
    batchSize: c_int,
) -> cusolverStatus_t
Expand description

The helper functions below can calculate the sizes needed for pre-allocated buffer.

The S and D data types are real valued single and double precision, respectively.

The C and Z data types are complex valued single and double precision, respectively.

This function gesvda (a stands for approximate) approximates the singular value decomposition of a tall skinny $m \times n$ matrix A and corresponding the left and right singular vectors. The economy form of SVD is written by: $$ A = U\*\Sigma\*V^{H} $$

where $\Sigma$ is an $n \times n$ matrix. U is an $m \times n$ unitary matrix, and V is an $n \times n$ unitary matrix. The diagonal elements of $\Sigma$ are the singular values of A; they are real and non-negative, and are returned in descending order. U and V are the left and right singular vectors of A.

gesvda computes eigenvalues of A**T*A, or A**H*A (if A is complex), to approximate singular values and singular vectors. It generates matrices U and V and transforms the matrix A to the following form: $$ U^{H}\*A\*V = S + E $$

where S is diagonal and E depends on rounding errors. To certain conditions, U, V and S approximate singular values and singular vectors up to machine zero of single precision. In general, V is unitary, S is more accurate than U. If singular value is far from zero, then left singular vector U is accurate. In other words, the accuracy of singular values and left singular vectors depend on the distance between singular value and zero. Since the computation of A**T*A, or A**H*A can greatly amplify errors, it is recommended to use gesvda only with well-conditioned data.

The input parameter rank decides the number of singular values and singular vectors are computed in parameter S, U and V.

The output parameter h_RnrmF computes Frobenius norm of residual. To compute h_RnrmF, info != NULL is required. $$ A - U\*S\*V^{H} $$

if the parameter rank is equal n. Otherwise, h_RnrmF reports: $$ {\|}U\*S\*V^{H}{\|} - {\|S\|} $$

in Frobenius norm sense, that is, how far U is from unitary.

gesvdaStridedBatched performs gesvda on each matrix. It requires that all matrices are of the same size m,n and are packed in a contiguous way, $$ \begin{split}A = \begin{pmatrix} {A0} & {A1} & \cdots \\ \end{pmatrix}\end{split} $$

Each matrix is column-major with leading dimension lda, so the formula for random access is $A_{k}\operatorname{(i,j)} = {A\lbrack\ i\ +\ lda\*j\ +\ strideA\*k\rbrack}$. Similarly, the formula for random access of S is $S_{k}\operatorname{(j)} = {S\lbrack\ j\ +\ StrideS\*k\rbrack}$, the formula for random access of U is $U_{k}\operatorname{(i,j)} = {U\lbrack\ i\ +\ ldu\*j\ +\ strideU\*k\rbrack}$ and the formula for random access of V is $V_{k}\operatorname{(i,j)} = {V\lbrack\ i\ +\ ldv\*j\ +\ strideV\*k\rbrack}$.

The user has to provide working space which is pointed by input parameter work. The input parameter lwork is the size of the working space, and it is returned by gesvdaStridedBatched_bufferSize(). Please note that the size in bytes of the working space is equal to sizeof(<type>) * lwork.

The output parameter info is an integer array of size batchSize. If the function returns cusolverStatus_t::CUSOLVER_STATUS_INVALID_VALUE, the first element info\[0\] = -i (less than zero) indicates i-th parameter is wrong (not counting handle). Otherwise, if info\[i\] = min(m,n)+1, gesvdaStridedBatched did not converge on the i-th matrix. If 0 < info\[i\] < min(m,n)+1, gesvdaStridedBatched could not compute an SVD of the i-th matrix fully; the leading singular values Si\[k\], 0 <= k <= info\[i\]-1, and corresponding singular vectors may still be useful. In this case, if h_RnrmF is requested, h_RnrmF reports the residual as if rank was set to info\[i\]-1.

Note that the problem size is limited by the condition batchSize*stride{A/S/U/V}<=INT32_MAX primarily due to the current implementation constraints.

Please visit cuSOLVER Library Samples - gesvdaStridedBatched for a code example.

Remark 1: The routine returns V, not $V^{H}$. This is different from gesvd.

Remark 2: The routine only supports m >=n.

Remark 3: It is recommended to use an FP64 data type, that is DgesvdaStridedBatched or ZgesvdaStridedBatched.

Remark 4: If the user is confident on the accuracy of singular values and singular vectors, for example, certain conditions hold (required singular value is far from zero), then the performance can be improved by passing a null pointer to h_RnrmF, i.e. no computation of the residual norm.