Intel® oneAPI Math Kernel Library Developer Reference - C
Computes groups of matrix-vector product with general matrices.
void cblas_sgemv_batch_strided (const CBLAS_LAYOUT layout, const CBLAS_TRANSPOSE trans, const MKL_INT m, const MKL_INT n, const float alpha, const float *a, const MKL_INT lda, const MKL_INT stridea, const float *x, const MKL_INT incx, const MKL_INT stridex, const float beta, float *y, const MKL_INT incy, const MKL_INT stridey, const MKL_INT batch_size);
void cblas_dgemv_batch_strided (const CBLAS_LAYOUT layout, const CBLAS_TRANSPOSE trans, const MKL_INT m, const MKL_INT n, const double alpha, const double *a, const MKL_INT lda, const MKL_INT stridea, const double *x, const MKL_INT incx, const MKL_INT stridex, const double beta, double *y, const MKL_INT incy, const MKL_INT stridey, const MKL_INT batch_size);
void cblas_cgemv_batch_strided (const CBLAS_LAYOUT layout, const CBLAS_TRANSPOSE trans, const MKL_INT m, const MKL_INT n, const void alpha, const void *a, const MKL_INT lda, const MKL_INT stridea, const void *x, const MKL_INT incx, const MKL_INT stridex, const void beta, void *y, const MKL_INT incy, const MKL_INT stridey, const MKL_INT batch_size);
void cblas_zgemv_batch_strided (const CBLAS_LAYOUT layout, const CBLAS_TRANSPOSE trans, const MKL_INT m, const MKL_INT n, const void alpha, const void *a, const MKL_INT lda, const MKL_INT stridea, const void *x, const MKL_INT incx, const MKL_INT stridex, const void beta, void *y, const MKL_INT incy, const MKL_INT stridey, const MKL_INT batch_size);
The cblas_?gemv_batch_strided routines perform a series of matrix-vector product added to a scaled vector. They are similar to the cblas_?gemv routine counterparts, but the cblas_?gemv_batch_strided routines perform matrix-vector operations with groups of matrices and vectors.
All matrices a and vectors x and y have the same parameters (size, increments) and are stored at constant stridea, stridex, and stridey from each other. The operation is defined as
for i = 0 … batch_size – 1 A is a matrix at offset i * stridea in a X and Y are vectors at offset i * stridex and i * stridey in x and y Y = alpha * op(A) * X + beta * Y end for
Specifies whether two-dimensional array storage is row-major (CblasRowMajor) or column-major (CblasColMajor).
Specifies op(A) the transposition operation applied to the A matrices.
if trans = CblasNoTrans, then op(A) = A;
if trans = CblasTrans, then op(A) = A';
if trans = CblasConjTrans, then op(A) = conjg(A').
Number of rows of the matrices A. The value of m must be at least 0.
Number of columns of the matrices A. The value of n must be at least 0.
Specifies the scalar alpha.
Array holding all the input matrix A. Must be of size at least lda*k + stridea * (batch_size -1) where k is n if column major layout is used or m if row major layout is used.
Specifies the leading dimension of the matrixA. It must be positive and at least mif column major layout is used or at least n if row major layout is used.
Stride between two consecutive A matrices, must be at least 0 if column major layout is used or at least n if row major layout is used.
Array holding all the input vector x. Must be of size at least (1 + (len-1)*abs(incx)) + stridex * (batch_size - 1) where len is n if the A matrix is not transposed or m otherwise.
Stride between two consecutive elements of the x vectors.
Stride between two consecutive x vectors, must be at least 0.
Specifies the scalar beta.
Array holding all the input vectors y. Must be of size at least batch_size * stridey.
Stride between two consecutive elements of the y vectors.
Stride between two consecutive y vectors, must be at least (1 + (len-1)*abs(incy)) where len is m if the matrix A is non transpose or n otherwise.
Number of gemv computations to perform and a matrices, x and y vectors. Must be at least 0.
Array holding the batch_size updated vector y.