Visible to Intel only — GUID: GUID-D60DAAB1-688C-4C87-8209-1DC20D324FC4
Visible to Intel only — GUID: GUID-D60DAAB1-688C-4C87-8209-1DC20D324FC4
gesvda_batch (USM Strided Version)
Computes the truncated SVD factorizations for a batch of general matrices. This routine belongs to the oneapi::mkl::lapack namespace.
Description
The routine computes the truncated SVD decomposition (or a reduced rank approximation) for a batch of general matrices Ai, as:
Ai = Ui * Si * ViT
where Ui and Vi are orthogonal matrices, Si is a diagonal matrix with singular values on the diagonal. Singular values are nonnegative and listed in decreasing order. A truncated SVD of a given matrix will produce matrices with the specified number of columns, where the number of columns is defined by the user or determined at runtime with the help of the user-defined tolerance threshold.
An approximation of each matrix can be also obtained as a product of two low-rank matrices (e.g. low rank product):
Ai = Pi * Qi
where
Pi = Ui * Si, Qi = ViT
if the number of rows is greater or equal to the number of columns,
Pi = Ui, Qi = Si` *ViT
otherwise.
The routine provides 3 possible ways for computing truncated SVD:
Compute truncated SVD with the help of the input array rank where rank[i] specifies the number of singular values and vectors to be computed in parameter Ui, Vi and Si for each matrix Ai.
Compute truncated SVD using a tolerance threshold. While computing SVD, singular values that are less than the user-defined tolerance are treated as zero, and they are not computed but set to zero.
Compute truncated SVD using the effective rank. The effective rank is determined by treating as zero those singular values that are less than the user-defined tolerance threshold times the largest singular value.
The routine can be also used for computing singular values only.
API
Syntax
namespace oneapi::mkl::lapack {
void gesvda_batch(sycl::queue &queue,
int64_t *iparm,
int64_t *irank,
int64_t m,
int64_t n,
T *a,
int64_t lda,
int64_t stride_a,
RealT *s,
int64_t stride_s,
T * u,
int64_t ldu,
int64_t stride_u,
T * vt,
int64_t ldvt,
int64_t stride_vt,
RealT tolerance,
RealT *residual,
int64_t batch_size,
T *scratchpad,
int64_t scratchpad_size,
const std::vector<sycl::event> &events = {})
}
gesvda_batch supports the following precisions and devices:
T |
Devices supported |
---|---|
float |
CPU, GPU* |
double |
CPU, GPU* |
std::complex<float> |
CPU, GPU* |
std::complex<double> |
CPU, GPU* |
*Interface support only; all computations are performed on the CPU.
Input Parameters
- queue
-
Device queue where calculations will be performed.
- iparm
-
Array of dimension 16 specifying options to compute truncated SVD. It also specifies the type of returned SVD decomposition form. The table below describes all individual components of the iparm parameter. Default values are denoted with an asterisk (*). NOTE: iparm[4]-iparm[15] are reserved for future use.
Component description |
Values |
---|---|
iparm[0] specifies a criterion for treating singular values as zeros |
If iparm[0]=-1, default iparm values (iparm[0-2]=0, iparm[3]=1) are used If iparm[0]=0*, the truncated SVD is computed with the help of the input array irank If iparm[0]=1, the truncated SVD is computed using the parameter tolerance If iparm[0]=2, the truncated SVD is computed using the effective rank. The effective rank is determined by treating as zero those singular values that are less than the user-defined tolerance times the largest singular value |
iparm[1] specifies option for computing singular vectors |
If iparm[1]=0*, both singular values and singular vectors are computed If iparm[1]=1, only singular values are computed |
iparm[2] specifies type of returned SVD decomposition |
If iparm[2]=0*, the truncated SVD as a product of three matrices Ai = Ui * Si * ViT is computed. If iparm[2]=1, the truncated as as a low rank product Ai = Pi * Qi is computed. |
iparm[3] specifies option for computing the residual vector |
If iparm[2]=0, the residual vector is not computed If iparm[2]=1*, the residual vector is computed |
- irank
-
If iparm[0]=0 or iparm[0]=-1, element irank[i] specifies the number of singular values and/or singular vectors to be computed in Ui, ViTand Si for each matrix Ai.
- m
-
The number of rows in the matrices Ai.
- n
-
The number of columns in the matrices Ai.
- a
-
Array holding input matrices Ai.
- lda
-
The leading dimension of a. Must be at least max(1, m).
- stride_a
-
The stride between the beginnings of matrices Ai inside the batch array a. Must be at least max(1, lda * n).
- stride_s
-
The stride between the beginnings of arrays Si inside the array s. Must be at least min(m,n).
- ldu
-
The leading dimension of Ui. Must be at least max(1, m).
- stride_u
-
The stride between the beginnings of matrices Ui inside the batch array u. Must be at least max(1, ldu * m).
- ldvt
-
The leading dimension of ViT . Must be at least max(1, n).
- stride_vt
-
The stride between the beginnings of matrices ViTinside the batch array vt. Must be at least max(1, ldvt * n).
- tolerance
-
Specifies the tolerance threshold. It is only used for computing truncated SVD in the cases of iparm[0]=1 and iparm[0]=2. Not used otherwise.
- batch_size
-
The number of problems in a batch.
- scratchpad
-
Scratchpad memory to be used by routine for storing intermediate results.
- scratchpad_size
-
Size of scratchpad memory as a number of floating point elements of type T. Size should not be less than the value returned by gesvda_batch_scratchpad_size (Strided Version) function.
Output Parameters
- irank
-
If iparm[0]=-1 or iparm[0]=0, element irank[i] is the number of computed singular values and/or singular vectors for matrix Ai.
- a
-
Unchanged on exit, if the residual vector is not required. Otherwise contains the residual matrix
Ai - Ui * Si * ViTif iparm[2]=0,
Ai - Pi * Qi if iparm[2]=1.
- s
-
Array of size min(m,n)*batch_size to store batch of singular values Si.
- u
-
u is array of size at least stride_u*batch_size to store batch of Ui if iparm[2]=0, or to store batch of Pi if iparm[2]=1.
- vt
-
vt is array of size at least stride_vt*batch_size to store batch of ViTif iparm[2]=0, or to store batch of Qi if iparm[2]=1.
- residual
-
Array of dimension batch_size. If iparm[3]=1, residual[i] is the Frobenius norm of the matrix:
|| Ai - Ui * Si * ViT ||
if iparm[2]=0, and
|| Ai - Pi * Qi ||
if iparm[2]=1.
Exceptions
Exception |
Description |
---|---|
mkl::lapack::batch_exception |
This exception is thrown when problems occur during calculations. You can obtain the info code of the problem using the info() method of the exception object: If info = -n, the n-th parameter had an illegal value. If info = 1, an internal memory allocation failed. If info = 2, an input parameter contains an invalid value. If info = 3, an error in algorithm while computing singular values . If info = 4, the routine encountered an empty structure or matrix array. |
Return Values
Output event to wait on to ensure computation is complete.