Developer Reference for Intel® oneAPI Math Kernel Library for C

ID 766684
Date 11/07/2023
Public

A newer version of this document is available. Customers should click here to go to the newest version.

Document Table of Contents

cblas_?gbmv

Computes a matrix-vector product with a general band matrix.

Syntax

void cblas_sgbmv (const CBLAS_LAYOUT Layout, const CBLAS_TRANSPOSE trans, const MKL_INT m, const MKL_INT n, const MKL_INT kl, const MKL_INT ku, const float alpha, const float *a, const MKL_INT lda, const float *x, const MKL_INT incx, const float beta, float *y, const MKL_INT incy);

void cblas_dgbmv (const CBLAS_LAYOUT Layout, const CBLAS_TRANSPOSE trans, const MKL_INT m, const MKL_INT n, const MKL_INT kl, const MKL_INT ku, const double alpha, const double *a, const MKL_INT lda, const double *x, const MKL_INT incx, const double beta, double *y, const MKL_INT incy);

void cblas_cgbmv (const CBLAS_LAYOUT Layout, const CBLAS_TRANSPOSE trans, const MKL_INT m, const MKL_INT n, const MKL_INT kl, const MKL_INT ku, const void *alpha, const void *a, const MKL_INT lda, const void *x, const MKL_INT incx, const void *beta, void *y, const MKL_INT incy);

void cblas_zgbmv (const CBLAS_LAYOUT Layout, const CBLAS_TRANSPOSE trans, const MKL_INT m, const MKL_INT n, const MKL_INT kl, const MKL_INT ku, const void *alpha, const void *a, const MKL_INT lda, const void *x, const MKL_INT incx, const void *beta, void *y, const MKL_INT incy);

Include Files

  • mkl.h

Description

The ?gbmv routines perform a matrix-vector operation defined as

y := alpha*A*x + beta*y,

or

y := alpha*A'*x + beta*y,

or

y := alpha *conjg(A')*x + beta*y,

where:

alpha and beta are scalars,

x and y are vectors,

A is an m-by-n band matrix, with kl sub-diagonals and ku super-diagonals.

Input Parameters

Layout

Specifies whether two-dimensional array storage is row-major (CblasRowMajor) or column-major (CblasColMajor).

trans

Specifies the operation:

If trans=CblasNoTrans, then y := alpha*A*x + beta*y

If trans=CblasTrans, then y := alpha*A'*x + beta*y

If trans=CblasConjTrans, then y := alpha *conjg(A')*x + beta*y

m

Specifies the number of rows of the matrix A.

The value of m must be at least zero.

n

Specifies the number of columns of the matrix A.

The value of n must be at least zero.

kl

Specifies the number of sub-diagonals of the matrix A.

The value of kl must satisfy 0kl.

ku

Specifies the number of super-diagonals of the matrix A.

The value of ku must satisfy 0ku.

alpha

Specifies the scalar alpha.

a

Array, size lda*n.

Layout = CblasColMajor: Before entry, the leading (kl + ku + 1) by n part of the array a must contain the matrix of coefficients. This matrix must be supplied column-by-column, with the leading diagonal of the matrix in row (ku) of the array, the first super-diagonal starting at position 1 in row (ku - 1), the first sub-diagonal starting at position 0 in row (ku + 1), and so on. Elements in the array a that do not correspond to elements in the band matrix (such as the top left ku by ku triangle) are not referenced.

The following program segment transfers a band matrix from conventional full matrix storage (matrix, with leading dimension ldm) to band storage (a, with leading dimension lda):

for (j = 0; j < n; j++) {
    k = ku - j;
    for (i = max(0, j-ku); i < min(m, j+kl+1); i++) {
        a[(k+i) + j*lda] = matrix[i + j*ldm];
    }
}

Layout = CblasRowMajor: Before entry, the leading (kl + ku + 1) by m part of the array a must contain the matrix of coefficients. This matrix must be supplied row-by-row, with the leading diagonal of the matrix in column (kl) of the array, the first super-diagonal starting at position 0 in column (kl + 1), the first sub-diagonal starting at position 1 in row (kl - 1), and so on. Elements in the array a that do not correspond to elements in the band matrix (such as the top left kl by kl triangle) are not referenced.

The following program segment transfers a band matrix from row-major full matrix storage (matrix, with leading dimension ldm) to band storage (a, with leading dimension lda):

for (i = 0; i < m; i++) {
    k = kl - i;
    for (j = max(0, i-kl); j < min(n, i+ku+1); j++) {
        a[(k+j) + i*lda] = matrix[j + i*ldm];
    }
}
lda

Specifies the leading dimension of a as declared in the calling (sub)program. The value of lda must be at least (kl + ku + 1).

x

Array, size at least (1 + (n - 1)*abs(incx)) when trans=CblasNoTrans, and at least (1 + (m - 1)*abs(incx)) otherwise. Before entry, the array x must contain the vector x.

incx

Specifies the increment for the elements of x. incx must not be zero.

beta

Specifies the scalar beta. When beta is equal to zero, then y need not be set on input.

y

Array, size at least (1 +(m - 1)*abs(incy)) when trans=CblasNoTrans and at least (1 +(n - 1)*abs(incy)) otherwise. Before entry, the incremented array y must contain the vector y.

incy

Specifies the increment for the elements of y.

The value of incy must not be zero.

Output Parameters

y

Buffer holding the updated vector y.