MAGMA  1.5.0
Matrix Algebra for GPU and Multicore Architectures
 All Functions Groups
double precision

Functions

void magmablas_cgemv (magma_trans_t trans, magma_int_t m, magma_int_t n, magmaFloatComplex alpha, const magmaFloatComplex *A, magma_int_t lda, const magmaFloatComplex *x, magma_int_t incx, magmaFloatComplex beta, magmaFloatComplex *y, magma_int_t incy)
 CGEMV performs one of the matrix-vector operations. More...
 
void magmablas_dgemv (magma_trans_t trans, magma_int_t m, magma_int_t n, double alpha, const double *A, magma_int_t lda, const double *x, magma_int_t incx, double beta, double *y, magma_int_t incy)
 This routine computes: 1) y = A x if trans == MagmaNoTrans, alpha == 1, beta == 0, and incx == incy == 1 (using magmablas code) 2) y = alpha A^T x if trans == MagmaTrans, beta == 0, and incx == incy == 1 (using magmablas code) 3) y = alpha A^trans x + beta y otherwise, using CUBLAS. More...
 
void magmablas_dgemv_MLU (magma_int_t m, magma_int_t n, const double *A, magma_int_t lda, const double *x, double *y)
 This routine computes y = y - Ax on the GPU. More...
 
void magmablas_dgemv_tesla (magma_trans_t trans, magma_int_t m, magma_int_t n, double alpha, const double *A, magma_int_t lda, const double *x, magma_int_t incx, double beta, double *y, magma_int_t incy)
 This routine computes: 1) y = A x if trans == MagmaNoTrans, alpha == 1, beta == 0, and incx == incy == 1 (using magmablas code) 2) y = alpha A^T x if trans == MagmaTrans, beta == 0, and incx == incy == 1 (using magmablas code) 3) y = alpha A^TRANS x + beta y otherwise, using CUBLAS. More...
 
magma_int_t magmablas_dsymv_work (magma_uplo_t uplo, magma_int_t n, double alpha, const double *A, magma_int_t lda, const double *x, magma_int_t incx, double beta, double *y, magma_int_t incy, double *dwork, magma_int_t lwork)
 magmablas_dsymv_work performs the matrix-vector operation: More...
 
magma_int_t magmablas_dsymv (magma_uplo_t uplo, magma_int_t n, double alpha, const double *A, magma_int_t lda, const double *x, magma_int_t incx, double beta, double *y, magma_int_t incy)
 magmablas_dsymv performs the matrix-vector operation: More...
 
magma_int_t magmablas_dsymv_mgpu_offset (magma_uplo_t uplo, magma_int_t n, double alpha, double **A, magma_int_t lda, double **x, magma_int_t incx, double beta, double **y, magma_int_t incy, double **work, magma_int_t lwork, magma_int_t num_gpus, magma_int_t nb, magma_int_t offset, magma_queue_t stream[][10])
 magmablas_dsymv performs the matrix-vector operation: More...
 
void magmablas_zgemv (magma_trans_t trans, magma_int_t m, magma_int_t n, magmaDoubleComplex alpha, const magmaDoubleComplex *A, magma_int_t lda, const magmaDoubleComplex *x, magma_int_t incx, magmaDoubleComplex beta, magmaDoubleComplex *y, magma_int_t incy)
 ZGEMV performs one of the matrix-vector operations. More...
 

Detailed Description

Function Documentation

void magmablas_cgemv ( magma_trans_t  trans,
magma_int_t  m,
magma_int_t  n,
magmaFloatComplex  alpha,
const magmaFloatComplex *  A,
magma_int_t  lda,
const magmaFloatComplex *  x,
magma_int_t  incx,
magmaFloatComplex  beta,
magmaFloatComplex *  y,
magma_int_t  incy 
)

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]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]ACOMPLEX array of dimension ( LDA, n ) on the GPU.
[in]ldaINTEGER LDA specifies the leading dimension of A.
[in]xCOMPLEX 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 REAL On entry, BETA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]yREAL 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.
void magmablas_dgemv ( magma_trans_t  trans,
magma_int_t  m,
magma_int_t  n,
double  alpha,
const double *  A,
magma_int_t  lda,
const double *  x,
magma_int_t  incx,
double  beta,
double *  y,
magma_int_t  incy 
)

This routine computes: 1) y = A x if trans == MagmaNoTrans, alpha == 1, beta == 0, and incx == incy == 1 (using magmablas code) 2) y = alpha A^T x if trans == MagmaTrans, beta == 0, and incx == incy == 1 (using magmablas code) 3) y = alpha A^trans x + beta y otherwise, using CUBLAS.

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^T*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 REAL On entry, ALPHA specifies the scalar alpha.
[in]ADOUBLE PRECISION array of dimension ( LDA, n ) on the GPU.
[in]ldaINTEGER LDA specifies the leading dimension of A.
[in]xDOUBLE 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 REAL On entry, BETA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]yDOUBLE 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.
void magmablas_dgemv_MLU ( magma_int_t  m,
magma_int_t  n,
const double *  A,
magma_int_t  lda,
const double *  x,
double *  y 
)

This routine computes y = y - Ax on the GPU.

Parameters
[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]ADOUBLE PRECISION array of dimension ( LDA, n ) on the GPU.
[in]ldaINTEGER. LDA specifies the leading dimension of A.
[in]xDOUBLE PRECISION array of dimension n.
[out]yDOUBLE PRECISION array of dimension n. On exit Y = Y - A X.
void magmablas_dgemv_tesla ( magma_trans_t  trans,
magma_int_t  m,
magma_int_t  n,
double  alpha,
const double *  A,
magma_int_t  lda,
const double *  x,
magma_int_t  incx,
double  beta,
double *  y,
magma_int_t  incy 
)

This routine computes: 1) y = A x if trans == MagmaNoTrans, alpha == 1, beta == 0, and incx == incy == 1 (using magmablas code) 2) y = alpha A^T x if trans == MagmaTrans, beta == 0, and incx == incy == 1 (using magmablas code) 3) y = alpha A^TRANS x + beta y otherwise, using CUBLAS.

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^T*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]ADOUBLE PRECISION array of dimension (LDA, N) on the GPU.
[in]ldaINTEGER LDA specifies the leading dimension of A.
[in]xDOUBLE PRECISION array of dimension n if trans == 'n' m if trans == 't'
[in]incxSpecifies the increment for the elements of X. INCX must not be zero.
[in]betaDOUBLE PRECISION On entry, BETA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]yDOUBLE PRECISION array of dimension m if trans == 'n' n if trans == 't'
[in]incySpecifies the increment for the elements of Y. INCY must not be zero.
magma_int_t magmablas_dsymv ( magma_uplo_t  uplo,
magma_int_t  n,
double  alpha,
const double *  A,
magma_int_t  lda,
const double *  x,
magma_int_t  incx,
double  beta,
double *  y,
magma_int_t  incy 
)

magmablas_dsymv performs the matrix-vector operation:

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

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

Parameters
[in]uplomagma_uplo_t. On entry, UPLO specifies whether the upper or lower triangular part of the array A is to be referenced as follows:
  • = MagmaUpper: Only the upper triangular part of A is to be referenced.
  • = MagmaLower: Only the lower triangular part of A is to be referenced.
[in]nINTEGER. On entry, N specifies the order of the matrix A. N must be at least zero.
[in]alphaDOUBLE PRECISION. On entry, ALPHA specifies the scalar alpha.
[in]ADOUBLE PRECISION array of DIMENSION ( LDA, n ). Before entry with UPLO = MagmaUpper, the leading n by n upper triangular part of the array A must contain the upper triangular part of the symmetric matrix and the strictly lower triangular part of A is not referenced. Before entry with UPLO = MagmaLower, the leading n by n lower triangular part of the array A must contain the lower triangular part of the symmetric matrix and the strictly upper triangular part of A is not referenced. Note that the imaginary parts of the diagonal elements need not be set and are assumed to be zero.
[in]ldaINTEGER. On entry, LDA specifies the first dimension of A as declared in the calling (sub) program. LDA must be at least max( 1, n ). It is recommended that lda is multiple of 16. Otherwise performance would be deteriorated as the memory accesses would not be fully coalescent.
[in]xDOUBLE PRECISION array of dimension at least ( 1 + ( n - 1 )*abs( INCX ) ). Before entry, the incremented array X must contain the n element vector x.
[in]incxINTEGER. On entry, INCX specifies the increment for the elements of X. INCX must not be zero.
[in]betaDOUBLE PRECISION. On entry, BETA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[in,out]yDOUBLE PRECISION array of dimension at least ( 1 + ( n - 1 )*abs( INCY ) ). Before entry, the incremented array Y must contain the n element vector y. On exit, Y is overwritten by the updated vector y.
[in]incyINTEGER. On entry, INCY specifies the increment for the elements of Y. INCY must not be zero.
magma_int_t magmablas_dsymv_mgpu_offset ( magma_uplo_t  uplo,
magma_int_t  n,
double  alpha,
double **  A,
magma_int_t  lda,
double **  x,
magma_int_t  incx,
double  beta,
double **  y,
magma_int_t  incy,
double **  work,
magma_int_t  lwork,
magma_int_t  num_gpus,
magma_int_t  nb,
magma_int_t  offset,
magma_queue_t  stream[][10] 
)

magmablas_dsymv performs the matrix-vector operation:

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

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

Parameters
[in]uplomagma_uplo_t. On entry, UPLO specifies whether the upper or lower triangular part of the array A is to be referenced as follows:
  • = MagmaUpper: Only the upper triangular part of A is to be referenced.
  • = MagmaLower: Only the lower triangular part of A is to be referenced.
[in]nINTEGER. On entry, N specifies the order of the matrix A. N must be at least zero.
[in]alphaDOUBLE PRECISION. On entry, ALPHA specifies the scalar alpha.
[in]ADOUBLE PRECISION array of DIMENSION ( LDA, n ). Before entry with UPLO = MagmaUpper, the leading n by n upper triangular part of the array A must contain the upper triangular part of the symmetric matrix and the strictly lower triangular part of A is not referenced. Before entry with UPLO = MagmaLower, the leading n by n lower triangular part of the array A must contain the lower triangular part of the symmetric matrix and the strictly upper triangular part of A is not referenced. Note that the imaginary parts of the diagonal elements need not be set and are assumed to be zero.
[in]ldaINTEGER. On entry, LDA specifies the first dimension of A as declared in the calling (sub) program. LDA must be at least max( 1, n ). It is recommended that lda is multiple of 16. Otherwise performance would be deteriorated as the memory accesses would not be fully coalescent.
[in]xDOUBLE PRECISION array of dimension at least ( 1 + ( n - 1 )*abs( INCX ) ). Before entry, the incremented array X must contain the n element vector x.
[in]incxINTEGER. On entry, INCX specifies the increment for the elements of X. INCX must not be zero.
[in]betaDOUBLE PRECISION. On entry, BETA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[in,out]yDOUBLE PRECISION array of dimension at least ( 1 + ( n - 1 )*abs( INCY ) ). Before entry, the incremented array Y must contain the n element vector y. On exit, Y is overwritten by the updated vector y.
[in]incyINTEGER. On entry, INCY specifies the increment for the elements of Y. INCY must not be zero.
magma_int_t magmablas_dsymv_work ( magma_uplo_t  uplo,
magma_int_t  n,
double  alpha,
const double *  A,
magma_int_t  lda,
const double *  x,
magma_int_t  incx,
double  beta,
double *  y,
magma_int_t  incy,
double *  dwork,
magma_int_t  lwork 
)

magmablas_dsymv_work performs the matrix-vector operation:

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

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

Parameters
[in]uplomagma_uplo_t. On entry, UPLO specifies whether the upper or lower triangular part of the array A is to be referenced as follows:
  • = MagmaUpper: Only the upper triangular part of A is to be referenced.
  • = MagmaLower: Only the lower triangular part of A is to be referenced.
[in]nINTEGER. On entry, N specifies the order of the matrix A. N must be at least zero.
[in]alphaDOUBLE PRECISION. On entry, ALPHA specifies the scalar alpha.
[in]ADOUBLE PRECISION array of DIMENSION ( LDA, n ). Before entry with UPLO = MagmaUpper, the leading n by n upper triangular part of the array A must contain the upper triangular part of the symmetric matrix and the strictly lower triangular part of A is not referenced. Before entry with UPLO = MagmaLower, the leading n by n lower triangular part of the array A must contain the lower triangular part of the symmetric matrix and the strictly upper triangular part of A is not referenced. Note that the imaginary parts of the diagonal elements need not be set and are assumed to be zero.
[in]ldaINTEGER. On entry, LDA specifies the first dimension of A as declared in the calling (sub) program. LDA must be at least max( 1, n ). It is recommended that lda is multiple of 16. Otherwise performance would be deteriorated as the memory accesses would not be fully coalescent.
[in]xDOUBLE PRECISION array of dimension at least ( 1 + ( n - 1 )*abs( INCX ) ). Before entry, the incremented array X must contain the n element vector x.
[in]incxINTEGER. On entry, INCX specifies the increment for the elements of X. INCX must not be zero.
[in]betaDOUBLE PRECISION. On entry, BETA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[in,out]yDOUBLE PRECISION array of dimension at least ( 1 + ( n - 1 )*abs( INCY ) ). Before entry, the incremented array Y must contain the n element vector y. On exit, Y is overwritten by the updated vector y.
[in]incyINTEGER. On entry, INCY specifies the increment for the elements of Y. INCY must not be zero.
[in]dwork(workspace) DOUBLE PRECISION array on the GPU, dimension (MAX(1, LWORK)),
[in]lworkINTEGER. The dimension of the array DWORK. LWORK >= LDA * ceil( N / NB_X ), where NB_X = 64.

MAGMA implements dsymv through two steps: 1) perform the multiplication in each thread block and put the intermediate value in dwork. 2) sum the intermediate values and store the final result in y.

magamblas_dsymv_work requires users to provide a workspace, while magmablas_dsymv is a wrapper routine allocating the workspace inside the routine and provides the same interface as cublas.

If users need to call dsymv frequently, we suggest using magmablas_dsymv_work instead of magmablas_dsymv. As the overhead to allocate and free in device memory in magmablas_dsymv would hurt performance. Our tests show that this penalty is about 10 Gflop/s when the matrix size is around 10000.

void magmablas_zgemv ( magma_trans_t  trans,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex  alpha,
const magmaDoubleComplex *  A,
magma_int_t  lda,
const magmaDoubleComplex *  x,
magma_int_t  incx,
magmaDoubleComplex  beta,
magmaDoubleComplex *  y,
magma_int_t  incy 
)

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]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]ACOMPLEX_16 array of dimension ( LDA, n ) on the GPU.
[in]ldaINTEGER LDA specifies the leading dimension of A.
[in]xCOMPLEX_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]betaDOUBLE REAL On entry, BETA specifies the scalar beta. When BETA is supplied as zero then Y need not be set on input.
[out]yDOUBLE 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.