MAGMA 2.9.0
Matrix Algebra for GPU and Multicore Architectures
Loading...
Searching...
No Matches

\(y = \alpha Ax + \beta y\) More...

Functions

void magmablas_cgemv_batched (magma_trans_t trans, magma_int_t m, magma_int_t n, const magmaFloatComplex alpha, magmaFloatComplex const *const *dA_array, magma_int_t ldda, magmaFloatComplex const *const *dx_array, magma_int_t incx, const magmaFloatComplex beta, magmaFloatComplex **dy_array, magma_int_t incy, magma_int_t batchCount, magma_queue_t queue)
 CGEMV performs one of the matrix-vector operations.
 
void magmablas_cgemv_batched_strided (magma_trans_t trans, magma_int_t m, magma_int_t n, const magmaFloatComplex alpha, const magmaFloatComplex *dA, magma_int_t ldda, magma_int_t strideA, const magmaFloatComplex *dx, magma_int_t incx, magma_int_t stridex, const magmaFloatComplex beta, magmaFloatComplex *dy, magma_int_t incy, magma_int_t stridey, magma_int_t batchCount, magma_queue_t queue)
 CGEMV performs one of the matrix-vector operations.
 
void magmablas_cgemv_vbatched (magma_trans_t trans, magma_int_t *m, magma_int_t *n, magmaFloatComplex alpha, magmaFloatComplex_ptr dA_array[], magma_int_t *ldda, magmaFloatComplex_ptr dx_array[], magma_int_t *incx, magmaFloatComplex beta, magmaFloatComplex_ptr dy_array[], magma_int_t *incy, magma_int_t batchCount, magma_queue_t queue)
 CGEMV performs one of the matrix-vector operations.
 
void magmablas_dgemv_batched (magma_trans_t trans, magma_int_t m, magma_int_t n, const double alpha, double const *const *dA_array, magma_int_t ldda, double const *const *dx_array, magma_int_t incx, const double beta, double **dy_array, magma_int_t incy, magma_int_t batchCount, magma_queue_t queue)
 DGEMV performs one of the matrix-vector operations.
 
void magmablas_dgemv_batched_strided (magma_trans_t trans, magma_int_t m, magma_int_t n, const double alpha, const double *dA, magma_int_t ldda, magma_int_t strideA, const double *dx, magma_int_t incx, magma_int_t stridex, const double beta, double *dy, magma_int_t incy, magma_int_t stridey, magma_int_t batchCount, magma_queue_t queue)
 DGEMV performs one of the matrix-vector operations.
 
void magmablas_dgemv_vbatched (magma_trans_t trans, magma_int_t *m, magma_int_t *n, double alpha, magmaDouble_ptr dA_array[], magma_int_t *ldda, magmaDouble_ptr dx_array[], magma_int_t *incx, double beta, magmaDouble_ptr dy_array[], magma_int_t *incy, magma_int_t batchCount, magma_queue_t queue)
 DGEMV performs one of the matrix-vector operations.
 
void magmablas_sgemv_batched (magma_trans_t trans, magma_int_t m, magma_int_t n, const float alpha, float const *const *dA_array, magma_int_t ldda, float const *const *dx_array, magma_int_t incx, const float beta, float **dy_array, magma_int_t incy, magma_int_t batchCount, magma_queue_t queue)
 SGEMV performs one of the matrix-vector operations.
 
void magmablas_sgemv_batched_strided (magma_trans_t trans, magma_int_t m, magma_int_t n, const float alpha, const float *dA, magma_int_t ldda, magma_int_t strideA, const float *dx, magma_int_t incx, magma_int_t stridex, const float beta, float *dy, magma_int_t incy, magma_int_t stridey, magma_int_t batchCount, magma_queue_t queue)
 SGEMV performs one of the matrix-vector operations.
 
void magmablas_sgemv_vbatched (magma_trans_t trans, magma_int_t *m, magma_int_t *n, float alpha, magmaFloat_ptr dA_array[], magma_int_t *ldda, magmaFloat_ptr dx_array[], magma_int_t *incx, float beta, magmaFloat_ptr dy_array[], magma_int_t *incy, magma_int_t batchCount, magma_queue_t queue)
 SGEMV performs one of the matrix-vector operations.
 
void magmablas_zgemv_batched (magma_trans_t trans, magma_int_t m, magma_int_t n, const magmaDoubleComplex alpha, magmaDoubleComplex const *const *dA_array, magma_int_t ldda, magmaDoubleComplex const *const *dx_array, magma_int_t incx, const magmaDoubleComplex beta, magmaDoubleComplex **dy_array, magma_int_t incy, magma_int_t batchCount, magma_queue_t queue)
 ZGEMV performs one of the matrix-vector operations.
 
void magmablas_zgemv_batched_strided (magma_trans_t trans, magma_int_t m, magma_int_t n, const magmaDoubleComplex alpha, const magmaDoubleComplex *dA, magma_int_t ldda, magma_int_t strideA, const magmaDoubleComplex *dx, magma_int_t incx, magma_int_t stridex, const magmaDoubleComplex beta, magmaDoubleComplex *dy, magma_int_t incy, magma_int_t stridey, magma_int_t batchCount, magma_queue_t queue)
 ZGEMV performs one of the matrix-vector operations.
 
void magmablas_zgemv_vbatched (magma_trans_t trans, magma_int_t *m, magma_int_t *n, magmaDoubleComplex alpha, magmaDoubleComplex_ptr dA_array[], magma_int_t *ldda, magmaDoubleComplex_ptr dx_array[], magma_int_t *incx, magmaDoubleComplex beta, magmaDoubleComplex_ptr dy_array[], magma_int_t *incy, magma_int_t batchCount, magma_queue_t queue)
 ZGEMV performs one of the matrix-vector operations.
 

Detailed Description

\(y = \alpha Ax + \beta y\)

Function Documentation

◆ magmablas_cgemv_batched()

void magmablas_cgemv_batched ( magma_trans_t trans,
magma_int_t m,
magma_int_t n,
const magmaFloatComplex alpha,
magmaFloatComplex const *const * dA_array,
magma_int_t ldda,
magmaFloatComplex const *const * dx_array,
magma_int_t incx,
const magmaFloatComplex beta,
magmaFloatComplex ** dy_array,
magma_int_t incy,
magma_int_t batchCount,
magma_queue_t queue )

CGEMV performs one of the matrix-vector operations.

y := alpha*A*x    + beta*y,   or
y := alpha*A**T*x + beta*y,   or
y := alpha*A**H*x + beta*y,

where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

This is the batch version of the routine, using pointer-to-pointer (P2P) interface. All matrices and vectors must have the same dimension(s).

Parameters
[in]transmagma_trans_t On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: y := alpha*A *x + beta*y
  • = MagmaTrans: y := alpha*A^T*x + beta*y
  • = MagmaConjTrans: y := alpha*A^H*x + beta*y
[in]mINTEGER On entry, m specifies the number of rows of the matrix A.
[in]nINTEGER On entry, n specifies the number of columns of the matrix A
[in]alphaCOMPLEX On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX array A of DIMENSION ( ldda, n ) on the GPU
[in]lddaINTEGER LDDA specifies the leading dimension of A.
[in]dx_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX array of dimension n if trans == MagmaNoTrans m if trans == MagmaTrans or MagmaConjTrans
[in]incxSpecifies the increment for the elements of X. INCX must not be zero.
[in]betaCOMPLEX On entry, ALPHA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]dy_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX array of dimension m if trans == MagmaNoTrans n if trans == MagmaTrans or MagmaConjTrans
[in]incySpecifies the increment for the elements of Y. INCY must not be zero.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.

◆ magmablas_cgemv_batched_strided()

void magmablas_cgemv_batched_strided ( magma_trans_t trans,
magma_int_t m,
magma_int_t n,
const magmaFloatComplex alpha,
const magmaFloatComplex * dA,
magma_int_t ldda,
magma_int_t strideA,
const magmaFloatComplex * dx,
magma_int_t incx,
magma_int_t stridex,
const magmaFloatComplex beta,
magmaFloatComplex * dy,
magma_int_t incy,
magma_int_t stridey,
magma_int_t batchCount,
magma_queue_t queue )

CGEMV performs one of the matrix-vector operations.

y := alpha*A*x    + beta*y,   or
y := alpha*A**T*x + beta*y,   or
y := alpha*A**H*x + beta*y,

where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

This is the batch version of the routine, using "pointer + stride" interface. All matrices and vectors must have the same dimension(s).

Parameters
[in]transmagma_trans_t On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: y := alpha*A *x + beta*y
  • = MagmaTrans: y := alpha*A^T*x + beta*y
  • = MagmaConjTrans: y := alpha*A^H*x + beta*y
[in]mINTEGER On entry, m specifies the number of rows of the matrix A.
[in]nINTEGER On entry, n specifies the number of columns of the matrix A
[in]alphaCOMPLEX On entry, ALPHA specifies the scalar alpha.
[in]dAPointer to the first COMPLEX array 'A' of DIMENSION ( ldda, n ) in the batch
[in]lddaINTEGER LDDA specifies the leading dimension of A.
[in]strideAINTEGER specifies the distance between two consecutive matrices in the batch.
[in]dxPointer to the first COMPLEX array 'x' in the batch, of dimension n if trans == MagmaNoTrans m if trans == MagmaTrans or MagmaConjTrans
[in]incxSpecifies the increment for the elements of X. INCX must not be zero.
[in]stridexINTEGER Specifies the distance between two consecutive vectors in the batch
[in]betaCOMPLEX On entry, ALPHA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]dyPointer to the first COMPLEX array 'y' in the batch, of dimension m if trans == MagmaNoTrans n if trans == MagmaTrans or MagmaConjTrans
[in]incySpecifies the increment for the elements of Y. INCY must not be zero.
[in]strideyINTEGER Specifies the distance between two consecutive vectors in the batch
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.

◆ magmablas_cgemv_vbatched()

void magmablas_cgemv_vbatched ( magma_trans_t trans,
magma_int_t * m,
magma_int_t * n,
magmaFloatComplex alpha,
magmaFloatComplex_ptr dA_array[],
magma_int_t * ldda,
magmaFloatComplex_ptr dx_array[],
magma_int_t * incx,
magmaFloatComplex beta,
magmaFloatComplex_ptr dy_array[],
magma_int_t * incy,
magma_int_t batchCount,
magma_queue_t queue )

CGEMV performs one of the matrix-vector operations.

y := alpha*A*x    + beta*y,   or
y := alpha*A**T*x + beta*y,   or
y := alpha*A**H*x + beta*y,

where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

Parameters
[in]transmagma_trans_t On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: y := alpha*A *x + beta*y
  • = MagmaTrans: y := alpha*A^T*x + beta*y
  • = MagmaConjTrans: y := alpha*A^H*x + beta*y
[in]mArray of integers, dimension (batchCount + 1). On entry, each INTEGER M specifies the number of rows of each matrix A. The last element of the array is used internally by the routine.
[in]nArray of integers, dimension (batchCount + 1). On entry, each INTEGER N specifies the number of columns of each matrix A The last element of the array is used internally by the routine.
[in]alphaCOMPLEX On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX array A of DIMENSION ( LDDA, N ) on the GPU
[in]lddaArray of integers, dimension (batchCount + 1). Each INTEGER LDDA specifies the leading dimension of each matrix A.
[in]dx_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX array of dimension N if trans == MagmaNoTrans M if trans == MagmaTrans or MagmaConjTrans
[in]incxArray of integers, dimension (batchCount + 1). Each integer specifies the increment for the elements of each vector X. INCX must not be zero. The last element of the array is used internally by the routine.
[in]betaCOMPLEX On entry, ALPHA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]dy_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX array of dimension M if trans == MagmaNoTrans N if trans == MagmaTrans or MagmaConjTrans
[in]incyArray of integers, dimension (batchCount + 1). Each integer specifies the increment for the elements of each vector Y. INCY must not be zero. The last element of the array is used internally by the routine.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.

◆ magmablas_dgemv_batched()

void magmablas_dgemv_batched ( magma_trans_t trans,
magma_int_t m,
magma_int_t n,
const double alpha,
double const *const * dA_array,
magma_int_t ldda,
double const *const * dx_array,
magma_int_t incx,
const double beta,
double ** dy_array,
magma_int_t incy,
magma_int_t batchCount,
magma_queue_t queue )

DGEMV performs one of the matrix-vector operations.

y := alpha*A*x    + beta*y,   or
y := alpha*A**T*x + beta*y,   or
y := alpha*A**H*x + beta*y,

where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

This is the batch version of the routine, using pointer-to-pointer (P2P) interface. All matrices and vectors must have the same dimension(s).

Parameters
[in]transmagma_trans_t On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: y := alpha*A *x + beta*y
  • = MagmaTrans: y := alpha*A^T*x + beta*y
  • = MagmaConjTrans: y := alpha*A^H*x + beta*y
[in]mINTEGER On entry, m specifies the number of rows of the matrix A.
[in]nINTEGER On entry, n specifies the number of columns of the matrix A
[in]alphaDOUBLE PRECISION On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayArray of pointers, dimension (batchCount). Each is a DOUBLE PRECISION array A of DIMENSION ( ldda, n ) on the GPU
[in]lddaINTEGER LDDA specifies the leading dimension of A.
[in]dx_arrayArray of pointers, dimension (batchCount). Each is a DOUBLE PRECISION array of dimension n if trans == MagmaNoTrans m if trans == MagmaTrans or MagmaConjTrans
[in]incxSpecifies the increment for the elements of X. INCX must not be zero.
[in]betaDOUBLE PRECISION On entry, ALPHA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]dy_arrayArray of pointers, dimension (batchCount). Each is a DOUBLE PRECISION array of dimension m if trans == MagmaNoTrans n if trans == MagmaTrans or MagmaConjTrans
[in]incySpecifies the increment for the elements of Y. INCY must not be zero.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.

◆ magmablas_dgemv_batched_strided()

void magmablas_dgemv_batched_strided ( magma_trans_t trans,
magma_int_t m,
magma_int_t n,
const double alpha,
const double * dA,
magma_int_t ldda,
magma_int_t strideA,
const double * dx,
magma_int_t incx,
magma_int_t stridex,
const double beta,
double * dy,
magma_int_t incy,
magma_int_t stridey,
magma_int_t batchCount,
magma_queue_t queue )

DGEMV performs one of the matrix-vector operations.

y := alpha*A*x    + beta*y,   or
y := alpha*A**T*x + beta*y,   or
y := alpha*A**H*x + beta*y,

where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

This is the batch version of the routine, using "pointer + stride" interface. All matrices and vectors must have the same dimension(s).

Parameters
[in]transmagma_trans_t On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: y := alpha*A *x + beta*y
  • = MagmaTrans: y := alpha*A^T*x + beta*y
  • = MagmaConjTrans: y := alpha*A^H*x + beta*y
[in]mINTEGER On entry, m specifies the number of rows of the matrix A.
[in]nINTEGER On entry, n specifies the number of columns of the matrix A
[in]alphaDOUBLE PRECISION On entry, ALPHA specifies the scalar alpha.
[in]dAPointer to the first DOUBLE PRECISION array 'A' of DIMENSION ( ldda, n ) in the batch
[in]lddaINTEGER LDDA specifies the leading dimension of A.
[in]strideAINTEGER specifies the distance between two consecutive matrices in the batch.
[in]dxPointer to the first DOUBLE PRECISION array 'x' in the batch, of dimension n if trans == MagmaNoTrans m if trans == MagmaTrans or MagmaConjTrans
[in]incxSpecifies the increment for the elements of X. INCX must not be zero.
[in]stridexINTEGER Specifies the distance between two consecutive vectors in the batch
[in]betaDOUBLE PRECISION On entry, ALPHA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]dyPointer to the first DOUBLE PRECISION array 'y' in the batch, of dimension m if trans == MagmaNoTrans n if trans == MagmaTrans or MagmaConjTrans
[in]incySpecifies the increment for the elements of Y. INCY must not be zero.
[in]strideyINTEGER Specifies the distance between two consecutive vectors in the batch
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.

◆ magmablas_dgemv_vbatched()

void magmablas_dgemv_vbatched ( magma_trans_t trans,
magma_int_t * m,
magma_int_t * n,
double alpha,
magmaDouble_ptr dA_array[],
magma_int_t * ldda,
magmaDouble_ptr dx_array[],
magma_int_t * incx,
double beta,
magmaDouble_ptr dy_array[],
magma_int_t * incy,
magma_int_t batchCount,
magma_queue_t queue )

DGEMV performs one of the matrix-vector operations.

y := alpha*A*x    + beta*y,   or
y := alpha*A**T*x + beta*y,   or
y := alpha*A**H*x + beta*y,

where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

Parameters
[in]transmagma_trans_t On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: y := alpha*A *x + beta*y
  • = MagmaTrans: y := alpha*A^T*x + beta*y
  • = MagmaConjTrans: y := alpha*A^H*x + beta*y
[in]mArray of integers, dimension (batchCount + 1). On entry, each INTEGER M specifies the number of rows of each matrix A. The last element of the array is used internally by the routine.
[in]nArray of integers, dimension (batchCount + 1). On entry, each INTEGER N specifies the number of columns of each matrix A The last element of the array is used internally by the routine.
[in]alphaDOUBLE PRECISION On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayArray of pointers, dimension (batchCount). Each is a DOUBLE PRECISION array A of DIMENSION ( LDDA, N ) on the GPU
[in]lddaArray of integers, dimension (batchCount + 1). Each INTEGER LDDA specifies the leading dimension of each matrix A.
[in]dx_arrayArray of pointers, dimension (batchCount). Each is a DOUBLE PRECISION array of dimension N if trans == MagmaNoTrans M if trans == MagmaTrans or MagmaConjTrans
[in]incxArray of integers, dimension (batchCount + 1). Each integer specifies the increment for the elements of each vector X. INCX must not be zero. The last element of the array is used internally by the routine.
[in]betaDOUBLE PRECISION On entry, ALPHA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]dy_arrayArray of pointers, dimension (batchCount). Each is a DOUBLE PRECISION array of dimension M if trans == MagmaNoTrans N if trans == MagmaTrans or MagmaConjTrans
[in]incyArray of integers, dimension (batchCount + 1). Each integer specifies the increment for the elements of each vector Y. INCY must not be zero. The last element of the array is used internally by the routine.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.

◆ magmablas_sgemv_batched()

void magmablas_sgemv_batched ( magma_trans_t trans,
magma_int_t m,
magma_int_t n,
const float alpha,
float const *const * dA_array,
magma_int_t ldda,
float const *const * dx_array,
magma_int_t incx,
const float beta,
float ** dy_array,
magma_int_t incy,
magma_int_t batchCount,
magma_queue_t queue )

SGEMV performs one of the matrix-vector operations.

y := alpha*A*x    + beta*y,   or
y := alpha*A**T*x + beta*y,   or
y := alpha*A**H*x + beta*y,

where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

This is the batch version of the routine, using pointer-to-pointer (P2P) interface. All matrices and vectors must have the same dimension(s).

Parameters
[in]transmagma_trans_t On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: y := alpha*A *x + beta*y
  • = MagmaTrans: y := alpha*A^T*x + beta*y
  • = MagmaConjTrans: y := alpha*A^H*x + beta*y
[in]mINTEGER On entry, m specifies the number of rows of the matrix A.
[in]nINTEGER On entry, n specifies the number of columns of the matrix A
[in]alphaREAL On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayArray of pointers, dimension (batchCount). Each is a REAL array A of DIMENSION ( ldda, n ) on the GPU
[in]lddaINTEGER LDDA specifies the leading dimension of A.
[in]dx_arrayArray of pointers, dimension (batchCount). Each is a REAL array of dimension n if trans == MagmaNoTrans m if trans == MagmaTrans or MagmaConjTrans
[in]incxSpecifies the increment for the elements of X. INCX must not be zero.
[in]betaREAL On entry, ALPHA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]dy_arrayArray of pointers, dimension (batchCount). Each is a REAL array of dimension m if trans == MagmaNoTrans n if trans == MagmaTrans or MagmaConjTrans
[in]incySpecifies the increment for the elements of Y. INCY must not be zero.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.

◆ magmablas_sgemv_batched_strided()

void magmablas_sgemv_batched_strided ( magma_trans_t trans,
magma_int_t m,
magma_int_t n,
const float alpha,
const float * dA,
magma_int_t ldda,
magma_int_t strideA,
const float * dx,
magma_int_t incx,
magma_int_t stridex,
const float beta,
float * dy,
magma_int_t incy,
magma_int_t stridey,
magma_int_t batchCount,
magma_queue_t queue )

SGEMV performs one of the matrix-vector operations.

y := alpha*A*x    + beta*y,   or
y := alpha*A**T*x + beta*y,   or
y := alpha*A**H*x + beta*y,

where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

This is the batch version of the routine, using "pointer + stride" interface. All matrices and vectors must have the same dimension(s).

Parameters
[in]transmagma_trans_t On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: y := alpha*A *x + beta*y
  • = MagmaTrans: y := alpha*A^T*x + beta*y
  • = MagmaConjTrans: y := alpha*A^H*x + beta*y
[in]mINTEGER On entry, m specifies the number of rows of the matrix A.
[in]nINTEGER On entry, n specifies the number of columns of the matrix A
[in]alphaREAL On entry, ALPHA specifies the scalar alpha.
[in]dAPointer to the first REAL array 'A' of DIMENSION ( ldda, n ) in the batch
[in]lddaINTEGER LDDA specifies the leading dimension of A.
[in]strideAINTEGER specifies the distance between two consecutive matrices in the batch.
[in]dxPointer to the first REAL array 'x' in the batch, of dimension n if trans == MagmaNoTrans m if trans == MagmaTrans or MagmaConjTrans
[in]incxSpecifies the increment for the elements of X. INCX must not be zero.
[in]stridexINTEGER Specifies the distance between two consecutive vectors in the batch
[in]betaREAL On entry, ALPHA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]dyPointer to the first REAL array 'y' in the batch, of dimension m if trans == MagmaNoTrans n if trans == MagmaTrans or MagmaConjTrans
[in]incySpecifies the increment for the elements of Y. INCY must not be zero.
[in]strideyINTEGER Specifies the distance between two consecutive vectors in the batch
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.

◆ magmablas_sgemv_vbatched()

void magmablas_sgemv_vbatched ( magma_trans_t trans,
magma_int_t * m,
magma_int_t * n,
float alpha,
magmaFloat_ptr dA_array[],
magma_int_t * ldda,
magmaFloat_ptr dx_array[],
magma_int_t * incx,
float beta,
magmaFloat_ptr dy_array[],
magma_int_t * incy,
magma_int_t batchCount,
magma_queue_t queue )

SGEMV performs one of the matrix-vector operations.

y := alpha*A*x    + beta*y,   or
y := alpha*A**T*x + beta*y,   or
y := alpha*A**H*x + beta*y,

where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

Parameters
[in]transmagma_trans_t On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: y := alpha*A *x + beta*y
  • = MagmaTrans: y := alpha*A^T*x + beta*y
  • = MagmaConjTrans: y := alpha*A^H*x + beta*y
[in]mArray of integers, dimension (batchCount + 1). On entry, each INTEGER M specifies the number of rows of each matrix A. The last element of the array is used internally by the routine.
[in]nArray of integers, dimension (batchCount + 1). On entry, each INTEGER N specifies the number of columns of each matrix A The last element of the array is used internally by the routine.
[in]alphaREAL On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayArray of pointers, dimension (batchCount). Each is a REAL array A of DIMENSION ( LDDA, N ) on the GPU
[in]lddaArray of integers, dimension (batchCount + 1). Each INTEGER LDDA specifies the leading dimension of each matrix A.
[in]dx_arrayArray of pointers, dimension (batchCount). Each is a REAL array of dimension N if trans == MagmaNoTrans M if trans == MagmaTrans or MagmaConjTrans
[in]incxArray of integers, dimension (batchCount + 1). Each integer specifies the increment for the elements of each vector X. INCX must not be zero. The last element of the array is used internally by the routine.
[in]betaREAL On entry, ALPHA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]dy_arrayArray of pointers, dimension (batchCount). Each is a REAL array of dimension M if trans == MagmaNoTrans N if trans == MagmaTrans or MagmaConjTrans
[in]incyArray of integers, dimension (batchCount + 1). Each integer specifies the increment for the elements of each vector Y. INCY must not be zero. The last element of the array is used internally by the routine.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.

◆ magmablas_zgemv_batched()

void magmablas_zgemv_batched ( magma_trans_t trans,
magma_int_t m,
magma_int_t n,
const magmaDoubleComplex alpha,
magmaDoubleComplex const *const * dA_array,
magma_int_t ldda,
magmaDoubleComplex const *const * dx_array,
magma_int_t incx,
const magmaDoubleComplex beta,
magmaDoubleComplex ** dy_array,
magma_int_t incy,
magma_int_t batchCount,
magma_queue_t queue )

ZGEMV performs one of the matrix-vector operations.

y := alpha*A*x    + beta*y,   or
y := alpha*A**T*x + beta*y,   or
y := alpha*A**H*x + beta*y,

where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

This is the batch version of the routine, using pointer-to-pointer (P2P) interface. All matrices and vectors must have the same dimension(s).

Parameters
[in]transmagma_trans_t On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: y := alpha*A *x + beta*y
  • = MagmaTrans: y := alpha*A^T*x + beta*y
  • = MagmaConjTrans: y := alpha*A^H*x + beta*y
[in]mINTEGER On entry, m specifies the number of rows of the matrix A.
[in]nINTEGER On entry, n specifies the number of columns of the matrix A
[in]alphaCOMPLEX_16 On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX_16 array A of DIMENSION ( ldda, n ) on the GPU
[in]lddaINTEGER LDDA specifies the leading dimension of A.
[in]dx_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX_16 array of dimension n if trans == MagmaNoTrans m if trans == MagmaTrans or MagmaConjTrans
[in]incxSpecifies the increment for the elements of X. INCX must not be zero.
[in]betaCOMPLEX_16 On entry, ALPHA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]dy_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX_16 array of dimension m if trans == MagmaNoTrans n if trans == MagmaTrans or MagmaConjTrans
[in]incySpecifies the increment for the elements of Y. INCY must not be zero.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.

◆ magmablas_zgemv_batched_strided()

void magmablas_zgemv_batched_strided ( magma_trans_t trans,
magma_int_t m,
magma_int_t n,
const magmaDoubleComplex alpha,
const magmaDoubleComplex * dA,
magma_int_t ldda,
magma_int_t strideA,
const magmaDoubleComplex * dx,
magma_int_t incx,
magma_int_t stridex,
const magmaDoubleComplex beta,
magmaDoubleComplex * dy,
magma_int_t incy,
magma_int_t stridey,
magma_int_t batchCount,
magma_queue_t queue )

ZGEMV performs one of the matrix-vector operations.

y := alpha*A*x    + beta*y,   or
y := alpha*A**T*x + beta*y,   or
y := alpha*A**H*x + beta*y,

where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

This is the batch version of the routine, using "pointer + stride" interface. All matrices and vectors must have the same dimension(s).

Parameters
[in]transmagma_trans_t On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: y := alpha*A *x + beta*y
  • = MagmaTrans: y := alpha*A^T*x + beta*y
  • = MagmaConjTrans: y := alpha*A^H*x + beta*y
[in]mINTEGER On entry, m specifies the number of rows of the matrix A.
[in]nINTEGER On entry, n specifies the number of columns of the matrix A
[in]alphaCOMPLEX_16 On entry, ALPHA specifies the scalar alpha.
[in]dAPointer to the first COMPLEX_16 array 'A' of DIMENSION ( ldda, n ) in the batch
[in]lddaINTEGER LDDA specifies the leading dimension of A.
[in]strideAINTEGER specifies the distance between two consecutive matrices in the batch.
[in]dxPointer to the first COMPLEX_16 array 'x' in the batch, of dimension n if trans == MagmaNoTrans m if trans == MagmaTrans or MagmaConjTrans
[in]incxSpecifies the increment for the elements of X. INCX must not be zero.
[in]stridexINTEGER Specifies the distance between two consecutive vectors in the batch
[in]betaCOMPLEX_16 On entry, ALPHA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]dyPointer to the first COMPLEX_16 array 'y' in the batch, of dimension m if trans == MagmaNoTrans n if trans == MagmaTrans or MagmaConjTrans
[in]incySpecifies the increment for the elements of Y. INCY must not be zero.
[in]strideyINTEGER Specifies the distance between two consecutive vectors in the batch
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.

◆ magmablas_zgemv_vbatched()

void magmablas_zgemv_vbatched ( magma_trans_t trans,
magma_int_t * m,
magma_int_t * n,
magmaDoubleComplex alpha,
magmaDoubleComplex_ptr dA_array[],
magma_int_t * ldda,
magmaDoubleComplex_ptr dx_array[],
magma_int_t * incx,
magmaDoubleComplex beta,
magmaDoubleComplex_ptr dy_array[],
magma_int_t * incy,
magma_int_t batchCount,
magma_queue_t queue )

ZGEMV performs one of the matrix-vector operations.

y := alpha*A*x    + beta*y,   or
y := alpha*A**T*x + beta*y,   or
y := alpha*A**H*x + beta*y,

where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

Parameters
[in]transmagma_trans_t On entry, TRANS specifies the operation to be performed as follows:
  • = MagmaNoTrans: y := alpha*A *x + beta*y
  • = MagmaTrans: y := alpha*A^T*x + beta*y
  • = MagmaConjTrans: y := alpha*A^H*x + beta*y
[in]mArray of integers, dimension (batchCount + 1). On entry, each INTEGER M specifies the number of rows of each matrix A. The last element of the array is used internally by the routine.
[in]nArray of integers, dimension (batchCount + 1). On entry, each INTEGER N specifies the number of columns of each matrix A The last element of the array is used internally by the routine.
[in]alphaCOMPLEX_16 On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX_16 array A of DIMENSION ( LDDA, N ) on the GPU
[in]lddaArray of integers, dimension (batchCount + 1). Each INTEGER LDDA specifies the leading dimension of each matrix A.
[in]dx_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX_16 array of dimension N if trans == MagmaNoTrans M if trans == MagmaTrans or MagmaConjTrans
[in]incxArray of integers, dimension (batchCount + 1). Each integer specifies the increment for the elements of each vector X. INCX must not be zero. The last element of the array is used internally by the routine.
[in]betaCOMPLEX_16 On entry, ALPHA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]dy_arrayArray of pointers, dimension (batchCount). Each is a COMPLEX_16 array of dimension M if trans == MagmaNoTrans N if trans == MagmaTrans or MagmaConjTrans
[in]incyArray of integers, dimension (batchCount + 1). Each integer specifies the increment for the elements of each vector Y. INCY must not be zero. The last element of the array is used internally by the routine.
[in]batchCountINTEGER The number of matrices to operate on.
[in]queuemagma_queue_t Queue to execute in.