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

Functions

magma_int_t magma_dcustomspmv (magma_int_t m, magma_int_t n, double alpha, double beta, double *x, double *y, magma_queue_t queue)
 This is an interface to any custom sparse matrix vector product. More...
 
magma_int_t magma_dge3pt (magma_int_t m, magma_int_t n, double alpha, double beta, magmaDouble_ptr dx, magmaDouble_ptr dy, magma_queue_t queue)
 This routine is a 3-pt-stencil operator derived from a FD-scheme in 2D with Dirichlet boundary. More...
 
magma_int_t magma_dgeaxpy (double alpha, magma_d_matrix X, double beta, magma_d_matrix *Y, magma_queue_t queue)
 This routine computes Y = alpha * X + beta * Y on the GPU. More...
 
magma_int_t magma_dgecsr5mv (magma_trans_t transA, magma_int_t m, magma_int_t n, magma_int_t p, double alpha, magma_int_t sigma, magma_int_t bit_y_offset, magma_int_t bit_scansum_offset, magma_int_t num_packet, magmaUIndex_ptr dtile_ptr, magmaUIndex_ptr dtile_desc, magmaIndex_ptr dtile_desc_offset_ptr, magmaIndex_ptr dtile_desc_offset, magmaDouble_ptr dcalibrator, magma_int_t tail_tile_start, magmaDouble_ptr dval, magmaIndex_ptr drowptr, magmaIndex_ptr dcolind, magmaDouble_ptr dx, double beta, magmaDouble_ptr dy, magma_queue_t queue)
 This routine computes y = alpha * A * x + beta * y on the GPU. More...
 
magma_int_t magma_dgecsrmv (magma_trans_t transA, magma_int_t m, magma_int_t n, double alpha, magmaDouble_ptr dval, magmaIndex_ptr drowptr, magmaIndex_ptr dcolind, magmaDouble_ptr dx, double beta, magmaDouble_ptr dy, magma_queue_t queue)
 This routine computes y = alpha * A * x + beta * y on the GPU. More...
 
magma_int_t magma_dgecsrmv_shift (magma_trans_t transA, magma_int_t m, magma_int_t n, double alpha, double lambda, magmaDouble_ptr dval, magmaIndex_ptr drowptr, magmaIndex_ptr dcolind, magmaDouble_ptr dx, double beta, magma_int_t offset, magma_int_t blocksize, magma_index_t *addrows, magmaDouble_ptr dy, magma_queue_t queue)
 This routine computes y = alpha * ( A -lambda I ) * x + beta * y on the GPU. More...
 
magma_int_t magma_dgecsrreimsplit (magma_d_matrix A, magma_d_matrix *ReA, magma_d_matrix *ImA, magma_queue_t queue)
 This routine takes an input matrix A in CSR format and located on the GPU and splits it into two matrixes ReA and ImA containing the real and the imaginary contributions of A. More...
 
magma_int_t magma_dgedensereimsplit (magma_d_matrix A, magma_d_matrix *ReA, magma_d_matrix *ImA, magma_queue_t queue)
 This routine takes an input matrix A in DENSE format and located on the GPU and splits it into two matrixes ReA and ImA containing the real and the imaginary contributions of A. More...
 
magma_int_t magma_dgeellmv (magma_trans_t transA, magma_int_t m, magma_int_t n, magma_int_t nnz_per_row, double alpha, magmaDouble_ptr dval, magmaIndex_ptr dcolind, magmaDouble_ptr dx, double beta, magmaDouble_ptr dy, magma_queue_t queue)
 This routine computes y = alpha * A * x + beta * y on the GPU. More...
 
magma_int_t magma_dgeellmv_shift (magma_trans_t transA, magma_int_t m, magma_int_t n, magma_int_t nnz_per_row, double alpha, double lambda, magmaDouble_ptr dval, magmaIndex_ptr dcolind, magmaDouble_ptr dx, double beta, magma_int_t offset, magma_int_t blocksize, magmaIndex_ptr addrows, magmaDouble_ptr dy, magma_queue_t queue)
 This routine computes y = alpha *( A - lambda I ) * x + beta * y on the GPU. More...
 
magma_int_t magma_dgeellrtmv (magma_trans_t transA, magma_int_t m, magma_int_t n, magma_int_t nnz_per_row, double alpha, magmaDouble_ptr dval, magmaIndex_ptr dcolind, magmaIndex_ptr drowlength, magmaDouble_ptr dx, double beta, magmaDouble_ptr dy, magma_int_t alignment, magma_int_t blocksize, magma_queue_t queue)
 This routine computes y = alpha * A * x + beta * y on the GPU. More...
 
magma_int_t magma_dgeelltmv_shift (magma_trans_t transA, magma_int_t m, magma_int_t n, magma_int_t nnz_per_row, double alpha, double lambda, magmaDouble_ptr dval, magmaIndex_ptr dcolind, magmaDouble_ptr dx, double beta, magma_int_t offset, magma_int_t blocksize, magmaIndex_ptr addrows, magmaDouble_ptr dy, magma_queue_t queue)
 This routine computes y = alpha *( A - lambda I ) * x + beta * y on the GPU. More...
 
magma_int_t magma_dmdotc (magma_int_t n, magma_int_t k, magmaDouble_ptr v, magmaDouble_ptr r, magmaDouble_ptr d1, magmaDouble_ptr d2, magmaDouble_ptr skp, magma_queue_t queue)
 Computes the scalar product of a set of vectors v_i such that. More...
 
magma_int_t magma_dgesellpmv (magma_trans_t transA, magma_int_t m, magma_int_t n, magma_int_t blocksize, magma_int_t slices, magma_int_t alignment, double alpha, magmaDouble_ptr dval, magmaIndex_ptr dcolind, magmaIndex_ptr drowptr, magmaDouble_ptr dx, double beta, magmaDouble_ptr dy, magma_queue_t queue)
 This routine computes y = alpha * A^t * x + beta * y on the GPU. More...
 
magma_int_t magma_dgesellcmv (magma_trans_t transA, magma_int_t m, magma_int_t n, magma_int_t blocksize, magma_int_t slices, magma_int_t alignment, double alpha, magmaDouble_ptr dval, magmaIndex_ptr dcolind, magmaIndex_ptr drowptr, magmaDouble_ptr dx, double beta, magmaDouble_ptr dy, magma_queue_t queue)
 This routine computes y = alpha * A^t * x + beta * y on the GPU. More...
 
magma_int_t magma_dmdotc_shfl (magma_int_t n, magma_int_t k, magmaDouble_ptr v, magmaDouble_ptr r, magmaDouble_ptr d1, magmaDouble_ptr d2, magmaDouble_ptr skp, magma_queue_t queue)
 Computes the scalar product of a set of vectors v_i such that. More...
 
magma_int_t magma_dmdotc4 (magma_int_t n, magmaDouble_ptr v0, magmaDouble_ptr w0, magmaDouble_ptr v1, magmaDouble_ptr w1, magmaDouble_ptr v2, magmaDouble_ptr w2, magmaDouble_ptr v3, magmaDouble_ptr w3, magmaDouble_ptr d1, magmaDouble_ptr d2, magmaDouble_ptr skp, magma_queue_t queue)
 Computes the scalar product of a set of 4 vectors such that. More...
 
magma_int_t magma_dmgecsrmv (magma_trans_t transA, magma_int_t m, magma_int_t n, magma_int_t num_vecs, double alpha, magmaDouble_ptr dval, magmaIndex_ptr drowptr, magmaIndex_ptr dcolind, magmaDouble_ptr dx, double beta, magmaDouble_ptr dy, magma_queue_t queue)
 This routine computes Y = alpha * A * X + beta * Y for X and Y sets of num_vec vectors on the GPU. More...
 
magma_int_t magma_dmgeellmv (magma_trans_t transA, magma_int_t m, magma_int_t n, magma_int_t num_vecs, magma_int_t nnz_per_row, double alpha, magmaDouble_ptr dval, magmaIndex_ptr dcolind, magmaDouble_ptr dx, double beta, magmaDouble_ptr dy, magma_queue_t queue)
 This routine computes Y = alpha * A * X + beta * Y for X and Y sets of num_vec vectors on the GPU. More...
 
magma_int_t magma_dmgeelltmv (magma_trans_t transA, magma_int_t m, magma_int_t n, magma_int_t num_vecs, magma_int_t nnz_per_row, double alpha, magmaDouble_ptr dval, magmaIndex_ptr dcolind, magmaDouble_ptr dx, double beta, magmaDouble_ptr dy, magma_queue_t queue)
 This routine computes Y = alpha * A * X + beta * Y for X and Y sets of num_vec vectors on the GPU. More...
 
magma_int_t magma_dmgesellpmv (magma_trans_t transA, magma_int_t m, magma_int_t n, magma_int_t num_vecs, magma_int_t blocksize, magma_int_t slices, magma_int_t alignment, double alpha, magmaDouble_ptr dval, magmaIndex_ptr dcolind, magmaIndex_ptr drowptr, magmaDouble_ptr dx, double beta, magmaDouble_ptr dy, magma_queue_t queue)
 This routine computes Y = alpha * A^t * X + beta * Y on the GPU. More...
 
magma_int_t magma_d_spmv (double alpha, magma_d_matrix A, magma_d_matrix x, double beta, magma_d_matrix y, magma_queue_t queue)
 For a given input matrix A and vectors x, y and scalars alpha, beta the wrapper determines the suitable SpMV computing y = alpha * A * x + beta * y. More...
 
magma_int_t magma_d_spmv_shift (double alpha, magma_d_matrix A, double lambda, magma_d_matrix x, double beta, magma_int_t offset, magma_int_t blocksize, magma_index_t *add_rows, magma_d_matrix y, magma_queue_t queue)
 For a given input matrix A and vectors x, y and scalars alpha, beta the wrapper determines the suitable SpMV computing y = alpha * ( A - lambda I ) * x + beta * y. More...
 
magma_int_t magma_d_spmm (double alpha, magma_d_matrix A, magma_d_matrix B, magma_d_matrix *C, magma_queue_t queue)
 For a given input matrix A and B and scalar alpha, the wrapper determines the suitable SpMV computing C = alpha * A * B. More...
 
magma_int_t magma_dcuspaxpy (double *alpha, magma_d_matrix A, double *beta, magma_d_matrix B, magma_d_matrix *AB, magma_queue_t queue)
 This is an interface to the cuSPARSE routine csrgeam computing the sum of two sparse matrices stored in csr format: More...
 
magma_int_t magma_dcuspmm (magma_d_matrix A, magma_d_matrix B, magma_d_matrix *AB, magma_queue_t queue)
 This is an interface to the cuSPARSE routine csrmm computing the product of two sparse matrices stored in csr format. More...
 

Detailed Description

Function Documentation

magma_int_t magma_dcustomspmv ( magma_int_t  m,
magma_int_t  n,
double  alpha,
double  beta,
double *  x,
double *  y,
magma_queue_t  queue 
)

This is an interface to any custom sparse matrix vector product.

It should compute y = alpha*FUNCTION(x) + beta*y The vectors are located on the device, the scalars on the CPU.

Parameters
[in]mmagma_int_t number of rows
[in]nmagma_int_t number of columns
[in]alphadouble scalar alpha
[in]xdouble * input vector x
[in]betadouble scalar beta
[out]ydouble * output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dge3pt ( magma_int_t  m,
magma_int_t  n,
double  alpha,
double  beta,
magmaDouble_ptr  dx,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine is a 3-pt-stencil operator derived from a FD-scheme in 2D with Dirichlet boundary.

It computes y_i = -2 x_i + x_{i-1} + x_{i+1}

Parameters
[in]mmagma_int_t number of rows in x and y
[in]nmagma_int_t number of columns in x and y
[in]alphadouble scalar multiplier
[in]betadouble scalar multiplier
[in]dxmagmaDouble_ptr input vector x
[out]dymagmaDouble_ptr output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dgeaxpy ( double  alpha,
magma_d_matrix  X,
double  beta,
magma_d_matrix *  Y,
magma_queue_t  queue 
)

This routine computes Y = alpha * X + beta * Y on the GPU.

The input format is magma_d_matrix. It can handle both, dense matrix (vector block) and CSR matrices. For the latter, it interfaces the cuSPARSE library.

Parameters
[in]alphadouble scalar multiplier.
[in]Xmagma_d_matrix input/output matrix Y.
[in]betadouble scalar multiplier.
[in,out]Ymagma_d_matrix* input matrix X.
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dgecsr5mv ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
magma_int_t  p,
double  alpha,
magma_int_t  sigma,
magma_int_t  bit_y_offset,
magma_int_t  bit_scansum_offset,
magma_int_t  num_packet,
magmaUIndex_ptr  dtile_ptr,
magmaUIndex_ptr  dtile_desc,
magmaIndex_ptr  dtile_desc_offset_ptr,
magmaIndex_ptr  dtile_desc_offset,
magmaDouble_ptr  dcalibrator,
magma_int_t  tail_tile_start,
magmaDouble_ptr  dval,
magmaIndex_ptr  drowptr,
magmaIndex_ptr  dcolind,
magmaDouble_ptr  dx,
double  beta,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine computes y = alpha * A * x + beta * y on the GPU.

The input format is CSR5 (val (tile-wise column-major), row_pointer, col (tile-wise column-major), tile_pointer, tile_desc).

Parameters
[in]transAmagma_trans_t transposition parameter for A
[in]mmagma_int_t number of rows in A
[in]nmagma_int_t number of columns in A
[in]pmagma_int_t number of tiles in A
[in]alphadouble scalar multiplier
[in]sigmamagma_int_t sigma in A in CSR5
[in]bit_y_offsetmagma_int_t bit_y_offset in A in CSR5
[in]bit_scansum_offsetmagma_int_t bit_scansum_offset in A in CSR5
[in]num_packetmagma_int_t num_packet in A in CSR5
[in]dtile_ptrmagmaUIndex_ptr tilepointer of A in CSR5
[in]dtile_descmagmaUIndex_ptr tiledescriptor of A in CSR5
[in]dtile_desc_offset_ptrmagmaIndex_ptr tiledescriptor_offsetpointer of A in CSR5
[in]dtile_desc_offsetmagmaIndex_ptr tiledescriptor_offsetpointer of A in CSR5
[in]dcalibratormagmaDouble_ptr calibrator of A in CSR5
[in]tail_tile_startmagma_int_t start of the last tile in A
[in]dvalmagmaDouble_ptr array containing values of A in CSR
[in]dvalmagmaDouble_ptr array containing values of A in CSR
[in]drowptrmagmaIndex_ptr rowpointer of A in CSR
[in]dcolindmagmaIndex_ptr columnindices of A in CSR
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar multiplier
[out]dymagmaDouble_ptr input/output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dgecsrmv ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
double  alpha,
magmaDouble_ptr  dval,
magmaIndex_ptr  drowptr,
magmaIndex_ptr  dcolind,
magmaDouble_ptr  dx,
double  beta,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine computes y = alpha * A * x + beta * y on the GPU.

The input format is CSR (val, row, col).

Parameters
[in]transAmagma_trans_t transposition parameter for A
[in]mmagma_int_t number of rows in A
[in]nmagma_int_t number of columns in A
[in]alphadouble scalar multiplier
[in]dvalmagmaDouble_ptr array containing values of A in CSR
[in]drowptrmagmaIndex_ptr rowpointer of A in CSR
[in]dcolindmagmaIndex_ptr columnindices of A in CSR
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar multiplier
[out]dymagmaDouble_ptr input/output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dgecsrmv_shift ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
double  alpha,
double  lambda,
magmaDouble_ptr  dval,
magmaIndex_ptr  drowptr,
magmaIndex_ptr  dcolind,
magmaDouble_ptr  dx,
double  beta,
magma_int_t  offset,
magma_int_t  blocksize,
magma_index_t *  addrows,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine computes y = alpha * ( A -lambda I ) * x + beta * y on the GPU.

It is a shifted version of the CSR-SpMV.

Parameters
[in]transAmagma_trans_t transposition parameter for A
[in]mmagma_int_t number of rows in A
[in]nmagma_int_t number of columns in A
[in]alphadouble scalar multiplier
[in]lambdadouble scalar multiplier
[in]dvalmagmaDouble_ptr array containing values of A in CSR
[in]drowptrmagmaIndex_ptr rowpointer of A in CSR
[in]dcolindmagmaIndex_ptr columnindices of A in CSR
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar multiplier
[in]offsetmagma_int_t in case not the main diagonal is scaled
[in]blocksizemagma_int_t in case of processing multiple vectors
[in]addrowsmagmaIndex_ptr in case the matrixpowerskernel is used
[out]dymagmaDouble_ptr output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dgecsrreimsplit ( magma_d_matrix  A,
magma_d_matrix *  ReA,
magma_d_matrix *  ImA,
magma_queue_t  queue 
)

This routine takes an input matrix A in CSR format and located on the GPU and splits it into two matrixes ReA and ImA containing the real and the imaginary contributions of A.

The output matrices are allocated within the routine.

Parameters
[in]Amagma_d_matrix input matrix A.
[out]ReAmagma_d_matrix* output matrix contaning real contributions.
[out]ImAmagma_d_matrix* output matrix contaning real contributions.
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dgedensereimsplit ( magma_d_matrix  A,
magma_d_matrix *  ReA,
magma_d_matrix *  ImA,
magma_queue_t  queue 
)

This routine takes an input matrix A in DENSE format and located on the GPU and splits it into two matrixes ReA and ImA containing the real and the imaginary contributions of A.

The output matrices are allocated within the routine.

Parameters
[in]Amagma_d_matrix input matrix A.
[out]ReAmagma_d_matrix* output matrix contaning real contributions.
[out]ImAmagma_d_matrix* output matrix contaning real contributions.
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dgeellmv ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
magma_int_t  nnz_per_row,
double  alpha,
magmaDouble_ptr  dval,
magmaIndex_ptr  dcolind,
magmaDouble_ptr  dx,
double  beta,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine computes y = alpha * A * x + beta * y on the GPU.

Input format is ELLPACK.

Parameters
[in]transAmagma_trans_t transposition parameter for A
[in]mmagma_int_t number of rows in A
[in]nmagma_int_t number of columns in A
[in]nnz_per_rowmagma_int_t number of elements in the longest row
[in]alphadouble scalar multiplier
[in]dvalmagmaDouble_ptr array containing values of A in ELLPACK
[in]dcolindmagmaIndex_ptr columnindices of A in ELLPACK
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar multiplier
[out]dymagmaDouble_ptr input/output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dgeellmv_shift ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
magma_int_t  nnz_per_row,
double  alpha,
double  lambda,
magmaDouble_ptr  dval,
magmaIndex_ptr  dcolind,
magmaDouble_ptr  dx,
double  beta,
magma_int_t  offset,
magma_int_t  blocksize,
magmaIndex_ptr  addrows,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine computes y = alpha *( A - lambda I ) * x + beta * y on the GPU.

Input format is ELLPACK. It is the shifted version of the ELLPACK SpMV.

Parameters
[in]transAmagma_trans_t transposition parameter for A
[in]mmagma_int_t number of rows in A
[in]nmagma_int_t number of columns in A
[in]nnz_per_rowmagma_int_t number of elements in the longest row
[in]alphadouble scalar multiplier
[in]lambdadouble scalar multiplier
[in]dvalmagmaDouble_ptr array containing values of A in ELLPACK
[in]dcolindmagmaIndex_ptr columnindices of A in ELLPACK
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar multiplier
[in]offsetmagma_int_t in case not the main diagonal is scaled
[in]blocksizemagma_int_t in case of processing multiple vectors
[in]addrowsmagmaIndex_ptr in case the matrixpowerskernel is used
[out]dymagmaDouble_ptr input/output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dgeellrtmv ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
magma_int_t  nnz_per_row,
double  alpha,
magmaDouble_ptr  dval,
magmaIndex_ptr  dcolind,
magmaIndex_ptr  drowlength,
magmaDouble_ptr  dx,
double  beta,
magmaDouble_ptr  dy,
magma_int_t  alignment,
magma_int_t  blocksize,
magma_queue_t  queue 
)

This routine computes y = alpha * A * x + beta * y on the GPU.

Input format is ELLRT. The ideas are taken from "Improving the performance of the sparse matrix vector product with GPUs", (CIT 2010), and modified to provide correct values.

Parameters
[in]transAmagma_trans_t transposition parameter for A
[in]mmagma_int_t number of rows
[in]nmagma_int_t number of columns
[in]nnz_per_rowmagma_int_t max number of nonzeros in a row
[in]alphadouble scalar alpha
[in]dvalmagmaDouble_ptr val array
[in]dcolindmagmaIndex_ptr col indices
[in]drowlengthmagmaIndex_ptr number of elements in each row
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar beta
[out]dymagmaDouble_ptr output vector y
[in]blocksizemagma_int_t threads per block
[in]alignmentmagma_int_t threads assigned to each row
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dgeelltmv_shift ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
magma_int_t  nnz_per_row,
double  alpha,
double  lambda,
magmaDouble_ptr  dval,
magmaIndex_ptr  dcolind,
magmaDouble_ptr  dx,
double  beta,
magma_int_t  offset,
magma_int_t  blocksize,
magmaIndex_ptr  addrows,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine computes y = alpha *( A - lambda I ) * x + beta * y on the GPU.

Input format is ELL.

Parameters
[in]transAmagma_trans_t transposition parameter for A
[in]mmagma_int_t number of rows in A
[in]nmagma_int_t number of columns in A
[in]nnz_per_rowmagma_int_t number of elements in the longest row
[in]alphadouble scalar multiplier
[in]lambdadouble scalar multiplier
[in]dvalmagmaDouble_ptr array containing values of A in ELL
[in]dcolindmagmaIndex_ptr columnindices of A in ELL
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar multiplier
[in]offsetmagma_int_t in case not the main diagonal is scaled
[in]blocksizemagma_int_t in case of processing multiple vectors
[in]addrowsmagmaIndex_ptr in case the matrixpowerskernel is used
[out]dymagmaDouble_ptr input/output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmdotc ( magma_int_t  n,
magma_int_t  k,
magmaDouble_ptr  v,
magmaDouble_ptr  r,
magmaDouble_ptr  d1,
magmaDouble_ptr  d2,
magmaDouble_ptr  skp,
magma_queue_t  queue 
)

Computes the scalar product of a set of vectors v_i such that.

skp = ( <v_0,r>, <v_1,r>, .. )

Returns the vector skp.

Parameters
[in]nint length of v_i and r
[in]kint

vectors v_i

Parameters
[in]vmagmaDouble_ptr v = (v_0 .. v_i.. v_k)
[in]rmagmaDouble_ptr r
[in]d1magmaDouble_ptr workspace
[in]d2magmaDouble_ptr workspace
[out]skpmagmaDouble_ptr vector[k] of scalar products (<v_i,r>...)
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dgesellpmv ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
magma_int_t  blocksize,
magma_int_t  slices,
magma_int_t  alignment,
double  alpha,
magmaDouble_ptr  dval,
magmaIndex_ptr  dcolind,
magmaIndex_ptr  drowptr,
magmaDouble_ptr  dx,
double  beta,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine computes y = alpha * A^t * x + beta * y on the GPU.

Input format is SELLP.

Parameters
[in]transAmagma_trans_t transposition parameter for A
[in]mmagma_int_t number of rows in A
[in]nmagma_int_t number of columns in A
[in]blocksizemagma_int_t number of rows in one ELL-slice
[in]slicesmagma_int_t number of slices in matrix
[in]alignmentmagma_int_t number of threads assigned to one row
[in]alphadouble scalar multiplier
[in]dvalmagmaDouble_ptr array containing values of A in SELLP
[in]dcolindmagmaIndex_ptr columnindices of A in SELLP
[in]drowptrmagmaIndex_ptr rowpointer of SELLP
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar multiplier
[out]dymagmaDouble_ptr input/output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dgesellcmv ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
magma_int_t  blocksize,
magma_int_t  slices,
magma_int_t  alignment,
double  alpha,
magmaDouble_ptr  dval,
magmaIndex_ptr  dcolind,
magmaIndex_ptr  drowptr,
magmaDouble_ptr  dx,
double  beta,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine computes y = alpha * A^t * x + beta * y on the GPU.

Input format is SELLC/SELLP.

Parameters
[in]transAmagma_trans_t transposition parameter for A
[in]mmagma_int_t number of rows in A
[in]nmagma_int_t number of columns in A
[in]blocksizemagma_int_t number of rows in one ELL-slice
[in]slicesmagma_int_t number of slices in matrix
[in]alignmentmagma_int_t number of threads assigned to one row (=1)
[in]alphadouble scalar multiplier
[in]dvalmagmaDouble_ptr array containing values of A in SELLC/P
[in]dcolindmagmaIndex_ptr columnindices of A in SELLC/P
[in]drowptrmagmaIndex_ptr rowpointer of SELLP
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar multiplier
[out]dymagmaDouble_ptr input/output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmdotc_shfl ( magma_int_t  n,
magma_int_t  k,
magmaDouble_ptr  v,
magmaDouble_ptr  r,
magmaDouble_ptr  d1,
magmaDouble_ptr  d2,
magmaDouble_ptr  skp,
magma_queue_t  queue 
)

Computes the scalar product of a set of vectors v_i such that.

skp = ( <v_0,r>, <v_1,r>, .. )

Returns the vector skp.

Parameters
[in]nint length of v_i and r
[in]kint

vectors v_i

Parameters
[in]vmagmaDouble_ptr v = (v_0 .. v_i.. v_k)
[in]rmagmaDouble_ptr r
[in]d1magmaDouble_ptr workspace
[in]d2magmaDouble_ptr workspace
[out]skpmagmaDouble_ptr vector[k] of scalar products (<v_i,r>...)
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmdotc4 ( magma_int_t  n,
magmaDouble_ptr  v0,
magmaDouble_ptr  w0,
magmaDouble_ptr  v1,
magmaDouble_ptr  w1,
magmaDouble_ptr  v2,
magmaDouble_ptr  w2,
magmaDouble_ptr  v3,
magmaDouble_ptr  w3,
magmaDouble_ptr  d1,
magmaDouble_ptr  d2,
magmaDouble_ptr  skp,
magma_queue_t  queue 
)

Computes the scalar product of a set of 4 vectors such that.

skp[0,1,2,3] = [ <v_0,w_0>, <v_1,w_1>, <v_2,w_2>, <v3,w_3> ]

Returns the vector skp. In case there are less dot products required, an easy workaround is given by doubling input.

Parameters
[in]nint length of v_i and w_i
[in]v0magmaDouble_ptr input vector
[in]w0magmaDouble_ptr input vector
[in]v1magmaDouble_ptr input vector
[in]w1magmaDouble_ptr input vector
[in]v2magmaDouble_ptr input vector
[in]w2magmaDouble_ptr input vector
[in]v3magmaDouble_ptr input vector
[in]w3magmaDouble_ptr input vector
[in]d1magmaDouble_ptr workspace
[in]d2magmaDouble_ptr workspace
[out]skpmagmaDouble_ptr vector[4] of scalar products [<v_i, w_i>] This vector is located on the host
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmgecsrmv ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
magma_int_t  num_vecs,
double  alpha,
magmaDouble_ptr  dval,
magmaIndex_ptr  drowptr,
magmaIndex_ptr  dcolind,
magmaDouble_ptr  dx,
double  beta,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine computes Y = alpha * A * X + beta * Y for X and Y sets of num_vec vectors on the GPU.

Input format is CSR.

Parameters
[in]transAmagma_trans_t transposition parameter for A
[in]mmagma_int_t number of rows in A
[in]nmagma_int_t number of columns in A
[in]num_vecsmama_int_t number of vectors
[in]alphadouble scalar multiplier
[in]dvalmagmaDouble_ptr array containing values of A in CSR
[in]drowptrmagmaIndex_ptr rowpointer of A in CSR
[in]dcolindmagmaIndex_ptr columnindices of A in CSR
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar multiplier
[out]dymagmaDouble_ptr input/output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmgeellmv ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
magma_int_t  num_vecs,
magma_int_t  nnz_per_row,
double  alpha,
magmaDouble_ptr  dval,
magmaIndex_ptr  dcolind,
magmaDouble_ptr  dx,
double  beta,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine computes Y = alpha * A * X + beta * Y for X and Y sets of num_vec vectors on the GPU.

Input format is ELLPACK.

Parameters
[in]transAmagma_trans_t transposition parameter for A
[in]mmagma_int_t number of rows in A
[in]nmagma_int_t number of columns in A
[in]num_vecsmama_int_t number of vectors
[in]nnz_per_rowmagma_int_t number of elements in the longest row
[in]alphadouble scalar multiplier
[in]dvalmagmaDouble_ptr array containing values of A in ELLPACK
[in]dcolindmagmaIndex_ptr columnindices of A in ELLPACK
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar multiplier
[out]dymagmaDouble_ptr input/output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmgeelltmv ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
magma_int_t  num_vecs,
magma_int_t  nnz_per_row,
double  alpha,
magmaDouble_ptr  dval,
magmaIndex_ptr  dcolind,
magmaDouble_ptr  dx,
double  beta,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine computes Y = alpha * A * X + beta * Y for X and Y sets of num_vec vectors on the GPU.

Input format is ELL.

Parameters
[in]transAmagma_trans_t transposition parameter for A
[in]mmagma_int_t number of rows in A
[in]nmagma_int_t number of columns in A
[in]num_vecsmama_int_t number of vectors
[in]nnz_per_rowmagma_int_t number of elements in the longest row
[in]alphadouble scalar multiplier
[in]dvalmagmaDouble_ptr array containing values of A in ELL
[in]dcolindmagmaIndex_ptr columnindices of A in ELL
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar multiplier
[out]dymagmaDouble_ptr input/output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dmgesellpmv ( magma_trans_t  transA,
magma_int_t  m,
magma_int_t  n,
magma_int_t  num_vecs,
magma_int_t  blocksize,
magma_int_t  slices,
magma_int_t  alignment,
double  alpha,
magmaDouble_ptr  dval,
magmaIndex_ptr  dcolind,
magmaIndex_ptr  drowptr,
magmaDouble_ptr  dx,
double  beta,
magmaDouble_ptr  dy,
magma_queue_t  queue 
)

This routine computes Y = alpha * A^t * X + beta * Y on the GPU.

Input format is SELLP. Note, that the input format for X is row-major while the output format for Y is column major!

Parameters
[in]transAmagma_trans_t transpose A?
[in]mmagma_int_t number of rows in A
[in]nmagma_int_t number of columns in A
[in]num_vecsmagma_int_t number of columns in X and Y
[in]blocksizemagma_int_t number of rows in one ELL-slice
[in]slicesmagma_int_t number of slices in matrix
[in]alignmentmagma_int_t number of threads assigned to one row
[in]alphadouble scalar multiplier
[in]dvalmagmaDouble_ptr array containing values of A in SELLP
[in]dcolindmagmaIndex_ptr columnindices of A in SELLP
[in]drowptrmagmaIndex_ptr rowpointer of SELLP
[in]dxmagmaDouble_ptr input vector x
[in]betadouble scalar multiplier
[out]dymagmaDouble_ptr input/output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_d_spmv ( double  alpha,
magma_d_matrix  A,
magma_d_matrix  x,
double  beta,
magma_d_matrix  y,
magma_queue_t  queue 
)

For a given input matrix A and vectors x, y and scalars alpha, beta the wrapper determines the suitable SpMV computing y = alpha * A * x + beta * y.

Parameters
[in]alphadouble scalar alpha
[in]Amagma_d_matrix sparse matrix A
[in]xmagma_d_matrix input vector x
[in]betadouble scalar beta
[out]ymagma_d_matrix output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_d_spmv_shift ( double  alpha,
magma_d_matrix  A,
double  lambda,
magma_d_matrix  x,
double  beta,
magma_int_t  offset,
magma_int_t  blocksize,
magma_index_t *  add_rows,
magma_d_matrix  y,
magma_queue_t  queue 
)

For a given input matrix A and vectors x, y and scalars alpha, beta the wrapper determines the suitable SpMV computing y = alpha * ( A - lambda I ) * x + beta * y.

Parameters
alphadouble scalar alpha
Amagma_d_matrix sparse matrix A
lambdadouble scalar lambda
xmagma_d_matrix input vector x
betadouble scalar beta
offsetmagma_int_t in case not the main diagonal is scaled
blocksizemagma_int_t in case of processing multiple vectors
add_rowsmagma_int_t* in case the matrixpowerskernel is used
ymagma_d_matrix output vector y
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_d_spmm ( double  alpha,
magma_d_matrix  A,
magma_d_matrix  B,
magma_d_matrix *  C,
magma_queue_t  queue 
)

For a given input matrix A and B and scalar alpha, the wrapper determines the suitable SpMV computing C = alpha * A * B.

Parameters
[in]alphadouble scalar alpha
[in]Amagma_d_matrix sparse matrix A
[in]Bmagma_d_matrix sparse matrix C
[out]Cmagma_d_matrix * outpur sparse matrix C
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcuspaxpy ( double *  alpha,
magma_d_matrix  A,
double *  beta,
magma_d_matrix  B,
magma_d_matrix *  AB,
magma_queue_t  queue 
)

This is an interface to the cuSPARSE routine csrgeam computing the sum of two sparse matrices stored in csr format:

C = alpha * A + beta * B
Parameters
[in]alphadouble* scalar
[in]Amagma_d_matrix input matrix
[in]betadouble* scalar
[in]Bmagma_d_matrix input matrix
[out]ABmagma_d_matrix* output matrix AB = alpha * A + beta * B
[in]queuemagma_queue_t Queue to execute in.
magma_int_t magma_dcuspmm ( magma_d_matrix  A,
magma_d_matrix  B,
magma_d_matrix *  AB,
magma_queue_t  queue 
)

This is an interface to the cuSPARSE routine csrmm computing the product of two sparse matrices stored in csr format.

Parameters
[in]Amagma_d_matrix input matrix
[in]Bmagma_d_matrix input matrix
[out]ABmagma_d_matrix* output matrix AB = A * B
[in]queuemagma_queue_t Queue to execute in.