Developer Reference for Intel® oneAPI Math Kernel Library for Fortran

ID 766686
Date 3/22/2024
Public
Document Table of Contents

?gesvda_batch_strided

Computes the truncated SVD of a group of general m-by-n matrices that are stored at a constant stride from each other in a contiguous block of memory.

Syntax

call sgesvda_batch_strided (iparm, irank, m, n, A, lda, stride_a, S, stride_s, U, ldu, stride_u, Vt, ldvt, stride_vt, tolerance, residual, work, lwork, batch_size, info )

call dgesvda_batch_strided (iparm, irank, m, n, A, lda, stride_a, S, stride_s, U, ldu, stride_u, Vt, ldvt, stride_vt, tolerance, residual, work, lwork, batch_size, info )

call cgesvda_batch_strided (iparm, irank, m, n, A, lda, stride_a, S, stride_s, U, ldu, stride_u, Vt, ldvt, stride_vt, tolerance, residual, work, lwork, batch_size, info )

call zgesvda_batch_strided (iparm, irank, m, n, A, lda, stride_a, S, stride_s, U, ldu, stride_u, Vt, ldvt, stride_vt, tolerance, residual, work, lwork, batch_size, info )

Include Files

mkl.fi

Description

The ?gesvda_batch_strided routines compute the truncated SVD for a group of general m-by-n matrices.

All matrices have the same parameters (matrix size, leading dimension) and are stored at constant stride_a from each other in a contiguous block of memory. The operation is defined as

for i = 0 … batch_size-1
    Ai is a matrix at offset i * stride_a from A
    Ai := Ui * Si*ViT
    Ai := Ui * Si *
end for

where Ui and Vi are orthogonal matrices, and 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 mxn matrix produces 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 (low-rank product):

Ai=Pi×Qi

where Pi=Ui×Si , Qi=ViT if m≥n, and Pi=Ui , Qi=Si × ViT otherwise.

The routines provide three possible ways to compute 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 parameters 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 of A is determined by treating as zero those singular values that are less than the user-defined tolerance threshold times the largest singular value.

The routines can be also used for computing singular values only.

Input Parameters

iparm

INTEGER. Array of dimension 16 specifying options to compute truncated SVD. Also specifies the type of returned SVD decomposition form. The individual components of the iparm parameter appear below. Default values are denoted with an asterisk (*).

iparm(1)

Specifies a criterion for treating singular values as zeros.

-1 Use default iparm values (iparm(1:3)=0, iparm(4)=1) . All other iparm settings are ignored.
= 0* Computes the truncated SVD with the help of the input array irank.
= 1 Computes the truncated SVD using the parameter tolerance.
= 2 Computes the truncated SVD using the effective rank. The effective rank of A is determined by treating as zero those singular values that are less than the user-defined tolerance multiplied by the largest singular value.

iparm(2)

Specifies the option for computing singular vectors.

0* Both singular values and singular vectors are computed.
1 Only singular values are computed.

iparm(3)

Specifies the type of the returned SVD decomposition.

0*

Computes the truncated SVD as a product of three matrices:

Ai=Ui×Si×ViT
1

Computes the truncated SVD as a low-rank product:

Ai=Pi×Qi

iparm(4)

Specifies the option for computing the residual vector.

0 The residual vector is not computed.
1* Computes the residual vector.
NOTE:

iparm(5)iparm(16) are reserved for future use.

irank

INTEGER. Array of size at least batch_size. If iparm(1)=0 or iparm(1)=-1, element irank(i) specifies the number of singular values and/or singular vectors to be computed in Ui , ViT, and Si for each matrix Ai.

m

INTEGER. The number of rows in the matrices Ai (m ≥ 0).

n

INTEGER. The number of columns in the matrices Ai (n ≥ 0).

A

REAL for sgesvda_batch_strided.

DOUBLE PRECISION for dgesvda_batch_strided.

COMPLEX for cgesvda_batch_strided.

DOUBLE COMPLEX for zgesvda_batch_strided.

Array of size at least stride_a * batch_size holding input matrices Ai.

lda

INTEGER. Specifies the leading dimension of the Ai matrices: lda ≥ max(1, m).

strde_a

INTEGER. Stride between two consecutive Ai matrices: stride_a ≥max(1, lda * n).

stride_s

INTEGER. The stride between two consecutive Si matrices: stride_s ≥ max(1, min(m,n)).

ldu

INTEGER. Specifies the leading dimension of the Ui matrices: ldu ≥ max(1, m).

stride_u

INTEGER. The stride between two consecutive Ui matrices: stride_u ≥ max(1, ldu * m).

ldvt

INTEGER. Specifies the leading dimension of the ViT matrices: ldvt ≥ max(1, n).

stride_vt

INTEGER. The stride between two consecutive ViT matrices: stride_vt ≥ max(1, ldvt * n).

tolerance

REAL for single precision flavors.

DOUBLE PRECISION for double precision flavors.

Specifies the tolerance threshold for computing truncated SVD in the cases of iparm(1)=1 and iparm(1)=2. Not used otherwise.

batch_size

INTEGER. The number of problems in a batch. Must be at least 0.

work

REAL for sgesvda_batch_strided.

DOUBLE PRECISION for dgesvda_batch_strided.

COMPLEX for cgesvda_batch_strided.

DOUBLE COMPLEX for zgesvda_batch_strided.

Workspace array with dimension max(1, lwork).

lwork

INTEGER. The dimension of the array work.

If lwork = -1, a workspace query is assumed: the routine only calculates the optimal size of the work array and returns this value as the first entry of the work array, and no error message related to lwork is issued by xerbla. If lwork is less than the required minimum size but is positive, the routine internally allocates the needed memory.

Output Parameters

irank

On exit, if iparm(1)=1 or iparm(1)=2, 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:=Ai- Ui×Si×ViT

if iparm(3)=0, and

A_i:=Ai-Pi×Qi

otherwise.

S

REAL for single precision flavors.

DOUBLE PRECISION for double precision flavors.

Array of size at least min(m,n)*batch_size to store a batch of singular values Si.

U

REAL for sgesvda_batch_strided.

DOUBLE PRECISION for dgesvda_batch_strided.

COMPLEX for cgesvda_batch_strided.

DOUBLE COMPLEX for zgesvda_batch_strided.

Array of size at least stride_u*batch_size to store a batch of Ui if iparm(3)=0, or to store a batch of Pi if iparm(3)=1.

Vt

REAL for sgesvda_batch_strided.

DOUBLE PRECISION for dgesvda_batch_strided.

COMPLEX for cgesvda_batch_strided.

DOUBLE COMPLEX for zgesvda_batch_strided.

Array of size at least stride_vt*batch_size to store a batch of ViT if iparm(3)=0, or to store a batch of Qi if iparm(3)=1.

residual

REAL for single precision flavors.

DOUBLE PRECISION for double precision flavors.

Array of dimension batch_size. If iparm(4)=1, residual(i) is the Frobenius norm of the matrix ||Ai - Ui×Si×ViT|| if iparm(3)=0, and ||Ai - Pi×Qi|| if iparm(3)=1.

info

INTEGER. Array of size at least batch_size, which reports the status for each matrix.

If info(i) = 0, the execution is successful for Ai.

If info(1) = -j, the j-th parameter had an illegal value.

If info(1) = 1, an internal memory allocation failed.

If info(i) = 2, an input parameter contains an invalid value.

If info(i) = 3, an error in algorithm while computing singular values of Ai occurred.

If info(1) = 4, the routine encountered an empty structure or matrix array.