MAGMA  1.6.1
Matrix Algebra for GPU and Multicore Architectures
 All Classes Files Functions Friends Groups Pages
double precision

Functions

void magma_dgemv (magma_trans_t transA, magma_int_t m, magma_int_t n, double alpha, magmaDouble_const_ptr dA, magma_int_t ldda, magmaDouble_const_ptr dx, magma_int_t incx, double beta, magmaDouble_ptr dy, magma_int_t incy)
 Perform matrix-vector product. More...
 
void magma_dger (magma_int_t m, magma_int_t n, double alpha, magmaDouble_const_ptr dx, magma_int_t incx, magmaDouble_const_ptr dy, magma_int_t incy, magmaDouble_ptr dA, magma_int_t ldda)
 Perform rank-1 update, \( A = \alpha x y^H + A \). More...
 
void magma_dsymv (magma_uplo_t uplo, magma_int_t n, double alpha, magmaDouble_const_ptr dA, magma_int_t ldda, magmaDouble_const_ptr dx, magma_int_t incx, double beta, magmaDouble_ptr dy, magma_int_t incy)
 Perform symmetric matrix-vector product, \( y = \alpha A x + \beta y \). More...
 
void magma_dsyr (magma_uplo_t uplo, magma_int_t n, double alpha, magmaDouble_const_ptr dx, magma_int_t incx, magmaDouble_ptr dA, magma_int_t ldda)
 Perform symmetric rank-1 update, \( A = \alpha x x^H + A \). More...
 
void magma_dsyr2 (magma_uplo_t uplo, magma_int_t n, double alpha, magmaDouble_const_ptr dx, magma_int_t incx, magmaDouble_const_ptr dy, magma_int_t incy, magmaDouble_ptr dA, magma_int_t ldda)
 Perform symmetric rank-2 update, \( A = \alpha x y^H + conj(\alpha) y x^H + A \). More...
 
void magma_dtrmv (magma_uplo_t uplo, magma_trans_t trans, magma_diag_t diag, magma_int_t n, magmaDouble_const_ptr dA, magma_int_t ldda, magmaDouble_ptr dx, magma_int_t incx)
 Perform triangular matrix-vector product. More...
 
void magma_dtrsv (magma_uplo_t uplo, magma_trans_t trans, magma_diag_t diag, magma_int_t n, magmaDouble_const_ptr dA, magma_int_t ldda, magmaDouble_ptr dx, magma_int_t incx)
 Solve triangular matrix-vector system (one right-hand side). More...
 
void magmablas_cgemv (magma_trans_t trans, magma_int_t m, magma_int_t n, magmaFloatComplex alpha, magmaFloatComplex_const_ptr dA, magma_int_t ldda, magmaFloatComplex_const_ptr dx, magma_int_t incx, magmaFloatComplex beta, magmaFloatComplex_ptr dy, magma_int_t incy)
 CGEMV performs one of the matrix-vector operations. More...
 
void magmablas_dgemv_batched (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)
 This routine computes Y = alpha opt(A) x + beta y, on the GPU, where A = dA_array[i],x = x_array[i] and y = y_array[i], i=[0,batchCount-1]. More...
 
void magmablas_dgemv_conjv (magma_int_t m, magma_int_t n, double alpha, magmaDouble_const_ptr dA, magma_int_t ldda, magmaDouble_const_ptr dx, magma_int_t incx, double beta, magmaDouble_ptr dy, magma_int_t incy)
 DGEMV_CONJV performs the matrix-vector operation. More...
 
void magmablas_dgemv (magma_trans_t trans, magma_int_t m, magma_int_t n, double alpha, magmaDouble_const_ptr dA, magma_int_t ldda, magmaDouble_const_ptr dx, magma_int_t incx, double beta, magmaDouble_ptr dy, magma_int_t incy)
 DGEMV performs one of the matrix-vector operations. 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...
 
void magmablas_dswapblk_q (magma_order_t order, magma_int_t n, magmaDouble_ptr dA, magma_int_t ldda, magmaDouble_ptr dB, magma_int_t lddb, magma_int_t i1, magma_int_t i2, const magma_int_t *ipiv, magma_int_t inci, magma_int_t offset, magma_queue_t queue)
 
 
 
void magmablas_dswapblk (magma_order_t order, magma_int_t n, magmaDouble_ptr dA, magma_int_t ldda, magmaDouble_ptr dB, magma_int_t lddb, magma_int_t i1, magma_int_t i2, const magma_int_t *ipiv, magma_int_t inci, magma_int_t offset)
 
magma_int_t magmablas_dsymv_work (magma_uplo_t uplo, magma_int_t n, double alpha, magmaDouble_const_ptr dA, magma_int_t ldda, magmaDouble_const_ptr dx, magma_int_t incx, double beta, magmaDouble_ptr dy, magma_int_t incy, magmaDouble_ptr dwork, magma_int_t lwork, magma_queue_t queue)
 magmablas_dsymv_work performs the matrix-vector operation: More...
 
magma_int_t magmablas_dsymv (magma_uplo_t uplo, magma_int_t n, double alpha, magmaDouble_const_ptr dA, magma_int_t ldda, magmaDouble_const_ptr dx, magma_int_t incx, double beta, magmaDouble_ptr dy, magma_int_t incy)
 magmablas_dsymv performs the matrix-vector operation: More...
 
magma_int_t magmablas_dsymv_mgpu (magma_uplo_t uplo, magma_int_t n, double alpha, magmaDouble_const_ptr const d_lA[], magma_int_t ldda, magma_int_t offset, double const *x, magma_int_t incx, double beta, double *y, magma_int_t incy, double *hwork, magma_int_t lhwork, magmaDouble_ptr dwork[], magma_int_t ldwork, magma_int_t ngpu, magma_int_t nb, magma_queue_t queues[])
 magmablas_dsymv_mgpu performs the matrix-vector operation: More...
 
magma_int_t magmablas_dsymv_mgpu_sync (magma_uplo_t uplo, magma_int_t n, double alpha, magmaDouble_const_ptr const d_lA[], magma_int_t ldda, magma_int_t offset, double const *x, magma_int_t incx, double beta, double *y, magma_int_t incy, double *hwork, magma_int_t lhwork, magmaDouble_ptr dwork[], magma_int_t ldwork, magma_int_t ngpu, magma_int_t nb, magma_queue_t queues[])
 Synchronizes and acculumates final dsymv result. More...
 
void magmablas_sgemv (magma_trans_t trans, magma_int_t m, magma_int_t n, float alpha, magmaFloat_const_ptr dA, magma_int_t ldda, magmaFloat_const_ptr dx, magma_int_t incx, float beta, magmaFloat_ptr dy, magma_int_t incy)
 SGEMV performs one of the matrix-vector operations. More...
 
void magmablas_zgemv (magma_trans_t trans, magma_int_t m, magma_int_t n, magmaDoubleComplex alpha, magmaDoubleComplex_const_ptr dA, magma_int_t ldda, magmaDoubleComplex_const_ptr dx, magma_int_t incx, magmaDoubleComplex beta, magmaDoubleComplex_ptr dy, magma_int_t incy)
 ZGEMV performs one of the matrix-vector operations. More...
 

Detailed Description

Function Documentation

void magma_dgemv ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
double  alpha,
magmaDouble_const_ptr  dA,
magma_int_t  ldda,
magmaDouble_const_ptr  dx,
magma_int_t  incx,
double  beta,
magmaDouble_ptr  dy,
magma_int_t  incy 
)

Perform matrix-vector product.

\( y = \alpha A x + \beta y \) (transA == MagmaNoTrans), or
\( y = \alpha A^T x + \beta y \) (transA == MagmaTrans), or
\( y = \alpha A^H x + \beta y \) (transA == MagmaConjTrans).

Parameters
[in]transAOperation to perform on A.
[in]mNumber of rows of A. m >= 0.
[in]nNumber of columns of A. n >= 0.
[in]alphaScalar \( \alpha \)
[in]dADOUBLE_PRECISION array of dimension (ldda,n), ldda >= max(1,m). The m-by-n matrix A, on GPU device.
[in]lddaLeading dimension of dA.
[in]dxDOUBLE_PRECISION array on GPU device. If transA == MagmaNoTrans, the n element vector x of dimension (1 + (n-1)*incx);
otherwise, the m element vector x of dimension (1 + (m-1)*incx).
[in]incxStride between consecutive elements of dx. incx != 0.
[in]betaScalar \( \beta \)
[in,out]dyDOUBLE_PRECISION array on GPU device. If transA == MagmaNoTrans, the m element vector y of dimension (1 + (m-1)*incy);
otherwise, the n element vector y of dimension (1 + (n-1)*incy).
[in]incyStride between consecutive elements of dy. incy != 0.
void magma_dger ( magma_int_t  m,
magma_int_t  n,
double  alpha,
magmaDouble_const_ptr  dx,
magma_int_t  incx,
magmaDouble_const_ptr  dy,
magma_int_t  incy,
magmaDouble_ptr  dA,
magma_int_t  ldda 
)

Perform rank-1 update, \( A = \alpha x y^H + A \).

Parameters
[in]mNumber of rows of A. m >= 0.
[in]nNumber of columns of A. n >= 0.
[in]alphaScalar \( \alpha \)
[in]dxDOUBLE_PRECISION array on GPU device. The m element vector x of dimension (1 + (m-1)*incx).
[in]incxStride between consecutive elements of dx. incx != 0.
[in]dyDOUBLE_PRECISION array on GPU device. The n element vector y of dimension (1 + (n-1)*incy).
[in]incyStride between consecutive elements of dy. incy != 0.
[in,out]dADOUBLE_PRECISION array on GPU device. The m-by-n matrix A of dimension (ldda,n), ldda >= max(1,m).
[in]lddaLeading dimension of dA.
void magma_dsymv ( magma_uplo_t  uplo,
magma_int_t  n,
double  alpha,
magmaDouble_const_ptr  dA,
magma_int_t  ldda,
magmaDouble_const_ptr  dx,
magma_int_t  incx,
double  beta,
magmaDouble_ptr  dy,
magma_int_t  incy 
)

Perform symmetric matrix-vector product, \( y = \alpha A x + \beta y \).

Parameters
[in]uploWhether the upper or lower triangle of A is referenced.
[in]nNumber of rows and columns of A. n >= 0.
[in]alphaScalar \( \alpha \)
[in]dADOUBLE_PRECISION array of dimension (ldda,n), ldda >= max(1,n). The n-by-n matrix A, on GPU device.
[in]lddaLeading dimension of dA.
[in]dxDOUBLE_PRECISION array on GPU device. The m element vector x of dimension (1 + (m-1)*incx).
[in]incxStride between consecutive elements of dx. incx != 0.
[in]betaScalar \( \beta \)
[in,out]dyDOUBLE_PRECISION array on GPU device. The n element vector y of dimension (1 + (n-1)*incy).
[in]incyStride between consecutive elements of dy. incy != 0.
void magma_dsyr ( magma_uplo_t  uplo,
magma_int_t  n,
double  alpha,
magmaDouble_const_ptr  dx,
magma_int_t  incx,
magmaDouble_ptr  dA,
magma_int_t  ldda 
)

Perform symmetric rank-1 update, \( A = \alpha x x^H + A \).

Parameters
[in]uploWhether the upper or lower triangle of A is referenced.
[in]nNumber of rows and columns of A. n >= 0.
[in]alphaScalar \( \alpha \)
[in]dxDOUBLE_PRECISION array on GPU device. The n element vector x of dimension (1 + (n-1)*incx).
[in]incxStride between consecutive elements of dx. incx != 0.
[in,out]dADOUBLE_PRECISION array of dimension (ldda,n), ldda >= max(1,n). The n-by-n matrix A, on GPU device.
[in]lddaLeading dimension of dA.
void magma_dsyr2 ( magma_uplo_t  uplo,
magma_int_t  n,
double  alpha,
magmaDouble_const_ptr  dx,
magma_int_t  incx,
magmaDouble_const_ptr  dy,
magma_int_t  incy,
magmaDouble_ptr  dA,
magma_int_t  ldda 
)

Perform symmetric rank-2 update, \( A = \alpha x y^H + conj(\alpha) y x^H + A \).

Parameters
[in]uploWhether the upper or lower triangle of A is referenced.
[in]nNumber of rows and columns of A. n >= 0.
[in]alphaScalar \( \alpha \)
[in]dxDOUBLE_PRECISION array on GPU device. The n element vector x of dimension (1 + (n-1)*incx).
[in]incxStride between consecutive elements of dx. incx != 0.
[in]dyDOUBLE_PRECISION array on GPU device. The n element vector y of dimension (1 + (n-1)*incy).
[in]incyStride between consecutive elements of dy. incy != 0.
[in,out]dADOUBLE_PRECISION array of dimension (ldda,n), ldda >= max(1,n). The n-by-n matrix A, on GPU device.
[in]lddaLeading dimension of dA.
void magma_dtrmv ( magma_uplo_t  uplo,
magma_trans_t  trans,
magma_diag_t  diag,
magma_int_t  n,
magmaDouble_const_ptr  dA,
magma_int_t  ldda,
magmaDouble_ptr  dx,
magma_int_t  incx 
)

Perform triangular matrix-vector product.

\( x = A x \) (trans == MagmaNoTrans), or
\( x = A^T x \) (trans == MagmaTrans), or
\( x = A^H x \) (trans == MagmaConjTrans).

Parameters
[in]uploWhether the upper or lower triangle of A is referenced.
[in]transOperation to perform on A.
[in]diagWhether the diagonal of A is assumed to be unit or non-unit.
[in]nNumber of rows and columns of A. n >= 0.
[in]dADOUBLE_PRECISION array of dimension (ldda,n), ldda >= max(1,n). The n-by-n matrix A, on GPU device.
[in]lddaLeading dimension of dA.
[in]dxDOUBLE_PRECISION array on GPU device. The n element vector x of dimension (1 + (n-1)*incx).
[in]incxStride between consecutive elements of dx. incx != 0.
void magma_dtrsv ( magma_uplo_t  uplo,
magma_trans_t  trans,
magma_diag_t  diag,
magma_int_t  n,
magmaDouble_const_ptr  dA,
magma_int_t  ldda,
magmaDouble_ptr  dx,
magma_int_t  incx 
)

Solve triangular matrix-vector system (one right-hand side).

\( A x = b \) (trans == MagmaNoTrans), or
\( A^T x = b \) (trans == MagmaTrans), or
\( A^H x = b \) (trans == MagmaConjTrans).

Parameters
[in]uploWhether the upper or lower triangle of A is referenced.
[in]transOperation to perform on A.
[in]diagWhether the diagonal of A is assumed to be unit or non-unit.
[in]nNumber of rows and columns of A. n >= 0.
[in]dADOUBLE_PRECISION array of dimension (ldda,n), ldda >= max(1,n). The n-by-n matrix A, on GPU device.
[in]lddaLeading dimension of dA.
[in,out]dxDOUBLE_PRECISION array on GPU device. On entry, the n element RHS vector b of dimension (1 + (n-1)*incx). On exit, overwritten with the solution vector x.
[in]incxStride between consecutive elements of dx. incx != 0.
void magmablas_cgemv ( magma_trans_t  trans,
magma_int_t  m,
magma_int_t  n,
magmaFloatComplex  alpha,
magmaFloatComplex_const_ptr  dA,
magma_int_t  ldda,
magmaFloatComplex_const_ptr  dx,
magma_int_t  incx,
magmaFloatComplex  beta,
magmaFloatComplex_ptr  dy,
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]dACOMPLEX array of dimension ( LDA, n ) on the GPU.
[in]ldaINTEGER LDA specifies the leading dimension of A.
[in]dxCOMPLEX 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]dyREAL 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,
magmaDouble_const_ptr  dA,
magma_int_t  ldda,
magmaDouble_const_ptr  dx,
magma_int_t  incx,
double  beta,
magmaDouble_ptr  dy,
magma_int_t  incy 
)

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]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]dADOUBLE_PRECISION array of dimension ( LDA, n ) on the GPU.
[in]ldaINTEGER LDA specifies the leading dimension of A.
[in]dxDOUBLE_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]dyDOUBLE 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_batched ( 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 
)

This routine computes Y = alpha opt(A) x + beta y, on the GPU, where A = dA_array[i],x = x_array[i] and y = y_array[i], i=[0,batchCount-1].

This is a batched version.

Parameters
[in]transCHARACTER*1. On entry, TRANS specifies the form of op( A ) to be used in the matrix multiplication as follows: = 'N': op( A ) = A. = 'T': op( A ) = A**T. = 'C': op( A ) = A**H.
[in]mINTEGER. On entry, M specifies the number of rows of the matrix opt(A).
[in]nINTEGER. On entry, N specifies the number of columns of the matrix opt(A)
[in]alphaDOUBLE PRECISION. On entry, ALPHA specifies the scalar alpha.
[in]dA_arrayA = dA_array[i] A: DOUBLE PRECISION array of dimension ( LDA, n ) on the GPU.
[in]ldaINTEGER. LDA specifies the leading dimension of A.
[in]x_arrayx = x_array[i] x: DOUBLE PRECISION array of dimension. n if trans == MagmaNoTrans. m if trans == MagmaTrans or MagmaConjTrans.
[in]incxINTEGER. incx specifies the increment for the elments of x. incx must not be zero.
[in]betaDOUBLE PRECISION. On entry, BETA specifies the scalar beta.
[out]y_arrayy = y_array[i]: On exit y = alpha opt(A) x + beta y. y: DOUBLE PRECISION array of dimension. m if trans == MagmaNoTrans. n if trans == MagmaTrans or MagmaConjTrans.
[in]incyINTEGER. incy specifies the increment for the elments of y. incy must not be zero.
[in]batchCountINTEGER number of pointers contained in dA_array, x_array and y_array.
void magmablas_dgemv_conjv ( magma_int_t  m,
magma_int_t  n,
double  alpha,
magmaDouble_const_ptr  dA,
magma_int_t  ldda,
magmaDouble_const_ptr  dx,
magma_int_t  incx,
double  beta,
magmaDouble_ptr  dy,
magma_int_t  incy 
)

DGEMV_CONJV performs the matrix-vector operation.

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

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

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]alphaDOUBLE_PRECISION On entry, ALPHA specifies the scalar alpha.
[in]dADOUBLE_PRECISION array of dimension ( LDA, n ) on the GPU.
[in]ldaINTEGER LDA specifies the leading dimension of A.
[in]dxDOUBLE_PRECISION array of dimension n
[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]dyDOUBLE PRECISION array of dimension m
[in]incySpecifies the increment for the elements of Y. INCY must not be zero.
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.
void magmablas_dswapblk ( magma_order_t  order,
magma_int_t  n,
magmaDouble_ptr  dA,
magma_int_t  ldda,
magmaDouble_ptr  dB,
magma_int_t  lddb,
magma_int_t  i1,
magma_int_t  i2,
const magma_int_t *  ipiv,
magma_int_t  inci,
magma_int_t  offset 
)
magma_int_t magmablas_dsymv ( magma_uplo_t  uplo,
magma_int_t  n,
double  alpha,
magmaDouble_const_ptr  dA,
magma_int_t  ldda,
magmaDouble_const_ptr  dx,
magma_int_t  incx,
double  beta,
magmaDouble_ptr  dy,
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]dADOUBLE_PRECISION array of DIMENSION ( LDDA, 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]lddaINTEGER. On entry, LDDA specifies the first dimension of A as declared in the calling (sub) program. LDDA must be at least max( 1, n ). It is recommended that ldda is multiple of 16. Otherwise performance would be deteriorated as the memory accesses would not be fully coalescent.
[in]dxDOUBLE_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]dyDOUBLE_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 ( magma_uplo_t  uplo,
magma_int_t  n,
double  alpha,
magmaDouble_const_ptr const  d_lA[],
magma_int_t  ldda,
magma_int_t  offset,
double const *  x,
magma_int_t  incx,
double  beta,
double *  y,
magma_int_t  incy,
double *  hwork,
magma_int_t  lhwork,
magmaDouble_ptr  dwork[],
magma_int_t  ldwork,
magma_int_t  ngpu,
magma_int_t  nb,
magma_queue_t  queues[] 
)

magmablas_dsymv_mgpu 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. Not currently supported.
  • = 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]d_lAArray of pointers, dimension (ngpu), to block-column distributed matrix A, with block size nb. d_lA[dev] is a DOUBLE_PRECISION array on GPU dev, of dimension (LDDA, nlocal), where
{ floor(n/nb/ngpu)*nb + nb if dev < floor(n/nb) % ngpu, nlocal = { floor(n/nb/ngpu)*nb + nnb if dev == floor(n/nb) % ngpu, { floor(n/nb/ngpu)*nb otherwise.
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]lddaINTEGER. On entry, LDDA specifies the first dimension of A as declared in the calling (sub) program. LDDA must be at least max( 1, n ). It is recommended that ldda is multiple of 16. Otherwise performance would be deteriorated as the memory accesses would not be fully coalescent.
[in]xDOUBLE_PRECISION array on the CPU (not the GPU), 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 on the CPU (not the GPU), 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.
hwork(workspace) DOUBLE_PRECISION array on the CPU, of dimension (lhwork).
[in]lhworkINTEGER. The dimension of the array hwork. lhwork >= ngpu*nb.
dwork(workspaces) Array of pointers, dimension (ngpu), to workspace on each GPU. dwork[dev] is a DOUBLE_PRECISION array on GPU dev, of dimension (ldwork).
[in]ldworkINTEGER. The dimension of each array dwork[dev]. ldwork >= ldda*( ceil((n + offset % nb) / nb) + 1 ).
[in]ngpuINTEGER. The number of GPUs to use.
[in]nbINTEGER. The block size used for distributing d_lA. Must be 64.
[in]queuesmagma_queue_t array of dimension (ngpu). queues[dev] is an execution queue on GPU dev.
magma_int_t magmablas_dsymv_mgpu_sync ( magma_uplo_t  uplo,
magma_int_t  n,
double  alpha,
magmaDouble_const_ptr const  d_lA[],
magma_int_t  ldda,
magma_int_t  offset,
double const *  x,
magma_int_t  incx,
double  beta,
double *  y,
magma_int_t  incy,
double *  hwork,
magma_int_t  lhwork,
magmaDouble_ptr  dwork[],
magma_int_t  ldwork,
magma_int_t  ngpu,
magma_int_t  nb,
magma_queue_t  queues[] 
)

Synchronizes and acculumates final dsymv result.

For convenience, the parameters are identical to magmablas_dsymv_mgpu (though some are unused here).

See also
magmablas_dsymv_mgpu
magma_int_t magmablas_dsymv_work ( magma_uplo_t  uplo,
magma_int_t  n,
double  alpha,
magmaDouble_const_ptr  dA,
magma_int_t  ldda,
magmaDouble_const_ptr  dx,
magma_int_t  incx,
double  beta,
magmaDouble_ptr  dy,
magma_int_t  incy,
magmaDouble_ptr  dwork,
magma_int_t  lwork,
magma_queue_t  queue 
)

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]dADOUBLE_PRECISION array of DIMENSION ( LDDA, 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]lddaINTEGER. On entry, LDDA specifies the first dimension of A as declared in the calling (sub) program. LDDA must be at least max( 1, n ). It is recommended that ldda is multiple of 16. Otherwise performance would be deteriorated as the memory accesses would not be fully coalescent.
[in]dxDOUBLE_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]dyDOUBLE_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 >= LDDA * ceil( N / NB_X ), where NB_X = 64.
[in]queuemagma_queue_t. Queue to execute in.

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_sgemv ( magma_trans_t  trans,
magma_int_t  m,
magma_int_t  n,
float  alpha,
magmaFloat_const_ptr  dA,
magma_int_t  ldda,
magmaFloat_const_ptr  dx,
magma_int_t  incx,
float  beta,
magmaFloat_ptr  dy,
magma_int_t  incy 
)

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]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]dAREAL array of dimension ( LDA, n ) on the GPU.
[in]ldaINTEGER LDA specifies the leading dimension of A.
[in]dxREAL 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]dyREAL 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_zgemv ( magma_trans_t  trans,
magma_int_t  m,
magma_int_t  n,
magmaDoubleComplex  alpha,
magmaDoubleComplex_const_ptr  dA,
magma_int_t  ldda,
magmaDoubleComplex_const_ptr  dx,
magma_int_t  incx,
magmaDoubleComplex  beta,
magmaDoubleComplex_ptr  dy,
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]dACOMPLEX_16 array of dimension ( LDA, n ) on the GPU.
[in]ldaINTEGER LDA specifies the leading dimension of A.
[in]dxCOMPLEX_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]dyDOUBLE 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.