Sparse Level 2 Functions#
This module holds all sparse level 2 routines.
The sparse level 2 routines describe operations between a matrix in sparse format and a vector in dense format.
rocsparse_bsrmv_ex_analysis()#
-
rocsparse_status rocsparse_sbsrmv_ex_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info)#
-
rocsparse_status rocsparse_dbsrmv_ex_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info)#
-
rocsparse_status rocsparse_cbsrmv_ex_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info)#
-
rocsparse_status rocsparse_zbsrmv_ex_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info)#
Sparse matrix vector multiplication using BSR storage format.
rocsparse_bsrmv_ex_analysis
performs the analysis step for rocsparse_sbsrmv(), rocsparse_dbsrmv(), rocsparse_cbsrmv() and rocsparse_zbsrmv(). It is expected that this function will be executed only once for a given matrix and particular operation type. The gathered analysis meta data can be cleared by rocsparse_bsrmv_ex_clear().Note
If the matrix sparsity pattern changes, the gathered information will become invalid.
Note
This function is blocking with respect to the host.
Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
dir – [in] matrix storage of BSR blocks.
trans – [in] matrix operation type.
mb – [in] number of block rows of the sparse BSR matrix.
nb – [in] number of block columns of the sparse BSR matrix.
nnzb – [in] number of non-zero blocks of the sparse BSR matrix.
descr – [in] descriptor of the sparse BSR matrix. Currently, only rocsparse_matrix_type_general is supported.
bsr_val – [in] array of
nnzb
blocks of the sparse BSR matrix.bsr_row_ptr – [in] array of
mb+1
elements that point to the start of every block row of the sparse BSR matrix.bsr_col_ind – [in] array of
nnzb
elements containing the block column indices of the sparse BSR matrix.block_dim – [in] block dimension of the sparse BSR matrix.
info – [out] structure that holds the information collected during the analysis step.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
mb
,nb
ornnzb
is invalid.rocsparse_status_invalid_pointer –
descr
,bsr_val
,bsr_row_ptr
,bsr_col_ind
orinfo
pointer is invalid.rocsparse_status_memory_error – the buffer for the gathered information could not be allocated.
rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_not_implemented –
trans
!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_bsrmv_ex()#
-
rocsparse_status rocsparse_sbsrmv_ex(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const float *alpha, const rocsparse_mat_descr descr, const float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, const float *x, const float *beta, float *y)#
-
rocsparse_status rocsparse_dbsrmv_ex(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const double *alpha, const rocsparse_mat_descr descr, const double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, const double *x, const double *beta, double *y)#
-
rocsparse_status rocsparse_cbsrmv_ex(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const rocsparse_float_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, const rocsparse_float_complex *x, const rocsparse_float_complex *beta, rocsparse_float_complex *y)#
-
rocsparse_status rocsparse_zbsrmv_ex(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const rocsparse_double_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, const rocsparse_double_complex *x, const rocsparse_double_complex *beta, rocsparse_double_complex *y)#
Sparse matrix vector multiplication using BSR storage format.
rocsparse_bsrmv_ex
multiplies the scalar \(\alpha\) with a sparse \((mb \cdot \text{block_dim}) \times (nb \cdot \text{block_dim})\) matrix, defined in BSR storage format, and the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == rocsparse_operation_none} \\ A^T, & \text{if trans == rocsparse_operation_transpose} \\ A^H, & \text{if trans == rocsparse_operation_conjugate_transpose} \end{array} \right. \end{split}\]Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
Currently, only
trans
== rocsparse_operation_none is supported.Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
dir – [in] matrix storage of BSR blocks.
trans – [in] matrix operation type.
mb – [in] number of block rows of the sparse BSR matrix.
nb – [in] number of block columns of the sparse BSR matrix.
nnzb – [in] number of non-zero blocks of the sparse BSR matrix.
alpha – [in] scalar \(\alpha\).
descr – [in] descriptor of the sparse BSR matrix. Currently, only rocsparse_matrix_type_general is supported.
bsr_val – [in] array of
nnzb
blocks of the sparse BSR matrix.bsr_row_ptr – [in] array of
mb+1
elements that point to the start of every block row of the sparse BSR matrix.bsr_col_ind – [in] array of
nnzb
elements containing the block column indices of the sparse BSR matrix.block_dim – [in] block dimension of the sparse BSR matrix.
x – [in] array of
nb*block_dim
elements ( \(op(A) = A\)) ormb*block_dim
elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).beta – [in] scalar \(\beta\).
y – [inout] array of
mb*block_dim
elements ( \(op(A) = A\)) ornb*block_dim
elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).info – [out] structure that holds the information collected during the analysis step.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
mb
,nb
,nnzb
orblock_dim
is invalid.rocsparse_status_invalid_pointer –
descr
,alpha
,bsr_val
,bsr_row_ind
,bsr_col_ind
,x
,beta
ory
pointer is invalid.rocsparse_status_arch_mismatch – the device is not supported.
rocsparse_status_not_implemented –
trans
!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_bsrmv_ex_clear()#
-
rocsparse_status rocsparse_bsrmv_ex_clear(rocsparse_handle handle, rocsparse_mat_info info)#
Sparse matrix vector multiplication using BSR storage format.
rocsparse_bsrmv_ex_clear
deallocates all memory that was allocated by rocsparse_sbsrmv_ex_analysis(), rocsparse_dbsrmv_ex_analysis(), rocsparse_cbsrmv_ex_analysis() or rocsparse_zbsrmv_ex_analysis(). This is especially useful, if memory is an issue and the analysis data is not required anymore for further computation, e.g. when switching to another sparse matrix format.Note
Calling
rocsparse_bsrmv_ex_clear
is optional. All allocated resources will be cleared, when the opaque rocsparse_mat_info struct is destroyed using rocsparse_destroy_mat_info().Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
info – [inout] structure that holds the information collected during analysis step.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_pointer –
info
pointer is invalid.rocsparse_status_memory_error – the buffer for the gathered information could not be deallocated.
rocsparse_status_internal_error – an internal error occurred.
rocsparse_bsrmv()#
-
rocsparse_status rocsparse_sbsrmv(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const float *alpha, const rocsparse_mat_descr descr, const float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, const float *x, const float *beta, float *y)#
-
rocsparse_status rocsparse_dbsrmv(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const double *alpha, const rocsparse_mat_descr descr, const double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, const double *x, const double *beta, double *y)#
-
rocsparse_status rocsparse_cbsrmv(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const rocsparse_float_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, const rocsparse_float_complex *x, const rocsparse_float_complex *beta, rocsparse_float_complex *y)#
-
rocsparse_status rocsparse_zbsrmv(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const rocsparse_double_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, const rocsparse_double_complex *x, const rocsparse_double_complex *beta, rocsparse_double_complex *y)#
Sparse matrix vector multiplication using BSR storage format.
rocsparse_bsrmv
multiplies the scalar \(\alpha\) with a sparse \((mb \cdot \text{block_dim}) \times (nb \cdot \text{block_dim})\) matrix, defined in BSR storage format, and the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == rocsparse_operation_none} \\ A^T, & \text{if trans == rocsparse_operation_transpose} \\ A^H, & \text{if trans == rocsparse_operation_conjugate_transpose} \end{array} \right. \end{split}\]- Example
This example performs a sparse matrix vector multiplication in BSR format.
// rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); // alpha * ( 1.0 0.0 2.0 ) * ( 1.0 ) + beta * ( 4.0 ) = ( 31.1 ) // ( 3.0 0.0 4.0 ) * ( 2.0 ) ( 5.0 ) = ( 62.0 ) // ( 5.0 6.0 0.0 ) * ( 3.0 ) ( 6.0 ) = ( 70.7 ) // ( 7.0 0.0 8.0 ) * ( 7.0 ) = ( 123.8 ) // BSR block dimension rocsparse_int bsr_dim = 2; // Number of block rows and columns rocsparse_int mb = 2; rocsparse_int nb = 2; // Number of non-zero blocks rocsparse_int nnzb = 4; // BSR row pointers rocsparse_int hbsr_row_ptr[3] = {0, 2, 4}; // BSR column indices rocsparse_int hbsr_col_ind[4] = {0, 1, 0, 1}; // BSR values double hbsr_val[16] = {1.0, 3.0, 0.0, 0.0, 2.0, 4.0, 0.0, 0.0, 5.0, 7.0, 6.0, 0.0, 0.0, 8.0, 0.0, 0.0}; // Block storage in column major rocsparse_direction dir = rocsparse_direction_column; // Transposition of the matrix rocsparse_operation trans = rocsparse_operation_none; // Scalar alpha and beta double alpha = 3.7; double beta = 1.3; // x and y double hx[4] = {1.0, 2.0, 3.0, 0.0}; double hy[4] = {4.0, 5.0, 6.0, 7.0}; // Matrix descriptor rocsparse_mat_descr descr; rocsparse_create_mat_descr(&descr); // Offload data to device rocsparse_int* dbsr_row_ptr; rocsparse_int* dbsr_col_ind; double* dbsr_val; double* dx; double* dy; hipMalloc((void**)&dbsr_row_ptr, sizeof(rocsparse_int) * (mb + 1)); hipMalloc((void**)&dbsr_col_ind, sizeof(rocsparse_int) * nnzb); hipMalloc((void**)&dbsr_val, sizeof(double) * nnzb * bsr_dim * bsr_dim); hipMalloc((void**)&dx, sizeof(double) * nb * bsr_dim); hipMalloc((void**)&dy, sizeof(double) * mb * bsr_dim); hipMemcpy(dbsr_row_ptr, hbsr_row_ptr, sizeof(rocsparse_int) * (mb + 1), hipMemcpyHostToDevice); hipMemcpy(dbsr_col_ind, hbsr_col_ind, sizeof(rocsparse_int) * nnzb, hipMemcpyHostToDevice); hipMemcpy(dbsr_val, hbsr_val, sizeof(double) * nnzb * bsr_dim * bsr_dim, hipMemcpyHostToDevice); hipMemcpy(dx, hx, sizeof(double) * nb * bsr_dim, hipMemcpyHostToDevice); hipMemcpy(dy, hy, sizeof(double) * mb * bsr_dim, hipMemcpyHostToDevice); rocsparse_mat_info info; rocsparse_create_mat_info(&info); // Call dbsrmv_analysis (Optional) rocsparse_dbsrmv_analysis(handle, dir, trans, mb, nb, nnzb, descr, dbsr_val, dbsr_row_ptr, dbsr_col_ind, bsr_dim, info); // Call dbsrmv to perform y = alpha * A x + beta * y rocsparse_dbsrmv(handle, dir, trans, mb, nb, nnzb, &alpha, descr, dbsr_val, dbsr_row_ptr, dbsr_col_ind, bsr_dim, info, dx, &beta, dy); // Copy result back to host hipMemcpy(hy, dy, sizeof(double) * mb * bsr_dim, hipMemcpyDeviceToHost); // Clear rocSPARSE rocsparse_destroy_mat_descr(descr); rocsparse_destroy_handle(handle); rocsparse_destroy_mat_info(info); // Clear device memory hipFree(dbsr_row_ptr); hipFree(dbsr_col_ind); hipFree(dbsr_val); hipFree(dx); hipFree(dy);
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
Currently, only
trans
== rocsparse_operation_none is supported.Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
dir – [in] matrix storage of BSR blocks.
trans – [in] matrix operation type.
mb – [in] number of block rows of the sparse BSR matrix.
nb – [in] number of block columns of the sparse BSR matrix.
nnzb – [in] number of non-zero blocks of the sparse BSR matrix.
alpha – [in] scalar \(\alpha\).
descr – [in] descriptor of the sparse BSR matrix. Currently, only rocsparse_matrix_type_general is supported.
bsr_val – [in] array of
nnzb
blocks of the sparse BSR matrix.bsr_row_ptr – [in] array of
mb+1
elements that point to the start of every block row of the sparse BSR matrix.bsr_col_ind – [in] array of
nnzb
elements containing the block column indices of the sparse BSR matrix.block_dim – [in] block dimension of the sparse BSR matrix.
x – [in] array of
nb*block_dim
elements ( \(op(A) = A\)) ormb*block_dim
elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).beta – [in] scalar \(\beta\).
y – [inout] array of
mb*block_dim
elements ( \(op(A) = A\)) ornb*block_dim
elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).info – [out] structure that holds the information collected during the analysis step.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
mb
,nb
,nnzb
orblock_dim
is invalid.rocsparse_status_invalid_pointer –
descr
,alpha
,bsr_val
,bsr_row_ind
,bsr_col_ind
,x
,beta
ory
pointer is invalid.rocsparse_status_arch_mismatch – the device is not supported.
rocsparse_status_not_implemented –
trans
!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_bsrxmv()#
-
rocsparse_status rocsparse_sbsrxmv(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int size_of_mask, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const float *alpha, const rocsparse_mat_descr descr, const float *bsr_val, const rocsparse_int *bsr_mask_ptr, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_end_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, const float *x, const float *beta, float *y)#
-
rocsparse_status rocsparse_dbsrxmv(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int size_of_mask, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const double *alpha, const rocsparse_mat_descr descr, const double *bsr_val, const rocsparse_int *bsr_mask_ptr, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_end_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, const double *x, const double *beta, double *y)#
-
rocsparse_status rocsparse_cbsrxmv(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int size_of_mask, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const rocsparse_float_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_mask_ptr, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_end_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, const rocsparse_float_complex *x, const rocsparse_float_complex *beta, rocsparse_float_complex *y)#
-
rocsparse_status rocsparse_zbsrxmv(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int size_of_mask, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const rocsparse_double_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_mask_ptr, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_end_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, const rocsparse_double_complex *x, const rocsparse_double_complex *beta, rocsparse_double_complex *y)#
Sparse matrix vector multiplication with mask operation using BSR storage format.
rocsparse_bsrxmv
multiplies the scalar \(\alpha\) with a sparse \((mb \cdot \text{block_dim}) \times (nb \cdot \text{block_dim})\) modified matrix, defined in BSR storage format, and the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \left( \alpha \cdot op(A) \cdot x + \beta \cdot y \right)\left( \text{mask} \right), \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == rocsparse_operation_none} \\ A^T, & \text{if trans == rocsparse_operation_transpose} \\ A^H, & \text{if trans == rocsparse_operation_conjugate_transpose} \end{array} \right. \end{split}\]The \(\text{mask}\) is defined as an array of block row indices. The input sparse matrix is defined with a modified BSR storage format where the beginning and the end of each row is defined with two arrays,
bsr_row_ptr
andbsr_end_ptr
(both of sizemb
), rather the usualbsr_row_ptr
of sizemb+1
.Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
Currently, only
trans
== rocsparse_operation_none is supported. Currently,block_dim==1
is not supported.Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
dir – [in] matrix storage of BSR blocks.
trans – [in] matrix operation type.
size_of_mask – [in] number of updated block rows of the array
y
.mb – [in] number of block rows of the sparse BSR matrix.
nb – [in] number of block columns of the sparse BSR matrix.
nnzb – [in] number of non-zero blocks of the sparse BSR matrix.
alpha – [in] scalar \(\alpha\).
descr – [in] descriptor of the sparse BSR matrix. Currently, only rocsparse_matrix_type_general is supported.
bsr_val – [in] array of
nnzb
blocks of the sparse BSR matrix.bsr_mask_ptr – [in] array of
size_of_mask
elements that give the indices of the updated block rows.bsr_row_ptr – [in] array of
mb
elements that point to the start of every block row of the sparse BSR matrix.bsr_end_ptr – [in] array of
mb
elements that point to the end of every block row of the sparse BSR matrix.bsr_col_ind – [in] array of
nnzb
elements containing the block column indices of the sparse BSR matrix.block_dim – [in] block dimension of the sparse BSR matrix.
x – [in] array of
nb*block_dim
elements ( \(op(A) = A\)) ormb*block_dim
elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).beta – [in] scalar \(\beta\).
y – [inout] array of
mb*block_dim
elements ( \(op(A) = A\)) ornb*block_dim
elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
mb
,nb
,nnzb
,block_dim
orsize_of_mask
is invalid.rocsparse_status_invalid_value –
size_of_mask
is greater thanmb
.rocsparse_status_invalid_pointer –
descr
,alpha
,bsr_val
,bsr_row_ind
,bsr_col_ind
,x
,beta
ory
pointer is invalid.rocsparse_status_arch_mismatch – the device is not supported.
rocsparse_status_not_implemented –
block_dim==1
,trans
!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_bsrsv_zero_pivot()#
-
rocsparse_status rocsparse_bsrsv_zero_pivot(rocsparse_handle handle, rocsparse_mat_info info, rocsparse_int *position)#
Sparse triangular solve using BSR storage format.
rocsparse_bsrsv_zero_pivot
returns rocsparse_status_zero_pivot, if either a structural or numerical zero has been found during rocsparse_sbsrsv_solve(), rocsparse_dbsrsv_solve(), rocsparse_cbsrsv_solve() or rocsparse_zbsrsv_solve() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored inposition
, using same index base as the BSR matrix.position
can be in host or device memory. If no zero pivot has been found,position
is set to -1 and rocsparse_status_success is returned instead.Note
rocsparse_bsrsv_zero_pivot
is a blocking function. It might influence performance negatively.Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
info – [in] structure that holds the information collected during the analysis step.
position – [inout] pointer to zero pivot \(j\), can be in host or device memory.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_pointer –
info
orposition
pointer is invalid.rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_zero_pivot – zero pivot has been found.
rocsparse_bsrsv_buffer_size()#
-
rocsparse_status rocsparse_sbsrsv_buffer_size(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, size_t *buffer_size)#
-
rocsparse_status rocsparse_dbsrsv_buffer_size(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, size_t *buffer_size)#
-
rocsparse_status rocsparse_cbsrsv_buffer_size(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, size_t *buffer_size)#
-
rocsparse_status rocsparse_zbsrsv_buffer_size(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, size_t *buffer_size)#
Sparse triangular solve using BSR storage format.
rocsparse_bsrsv_buffer_size
returns the size of the temporary storage buffer that is required by rocsparse_sbsrsv_analysis(), rocsparse_dbsrsv_analysis(), rocsparse_cbsrsv_analysis(), rocsparse_zbsrsv_analysis(), rocsparse_sbsrsv_solve(), rocsparse_dbsrsv_solve(), rocsparse_cbsrsv_solve() and rocsparse_zbsrsv_solve(). The temporary storage buffer must be allocated by the user.Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
dir – [in] matrix storage of BSR blocks.
trans – [in] matrix operation type.
mb – [in] number of block rows of the sparse BSR matrix.
nnzb – [in] number of non-zero blocks of the sparse BSR matrix.
descr – [in] descriptor of the sparse BSR matrix.
bsr_val – [in] array of
nnzb
blocks of the sparse BSR matrix.bsr_row_ptr – [in] array of
mb+1
elements that point to the start of every block row of the sparse BSR matrix.bsr_col_ind – [in] array of
nnz
containing the block column indices of the sparse BSR matrix.block_dim – [in] block dimension of the sparse BSR matrix.
info – [out] structure that holds the information collected during the analysis step.
buffer_size – [out] number of bytes of the temporary storage buffer required by rocsparse_sbsrsv_analysis(), rocsparse_dbsrsv_analysis(), rocsparse_cbsrsv_analysis(), rocsparse_zbsrsv_analysis(), rocsparse_sbsrsv_solve(), rocsparse_dbsrsv_solve(), rocsparse_cbsrsv_solve() and rocsparse_zbsrsv_solve().
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
mb
,nnzb
orblock_dim
is invalid.rocsparse_status_invalid_pointer –
descr
,bsr_val
,bsr_row_ptr
,bsr_col_ind
,info
orbuffer_size
pointer is invalid.rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_not_implemented –
trans
== rocsparse_operation_conjugate_transpose or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_bsrsv_analysis()#
-
rocsparse_status rocsparse_sbsrsv_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
-
rocsparse_status rocsparse_dbsrsv_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
-
rocsparse_status rocsparse_cbsrsv_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
-
rocsparse_status rocsparse_zbsrsv_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
Sparse triangular solve using BSR storage format.
rocsparse_bsrsv_analysis
performs the analysis step for rocsparse_sbsrsv_solve(), rocsparse_dbsrsv_solve(), rocsparse_cbsrsv_solve() and rocsparse_zbsrsv_solve(). It is expected that this function will be executed only once for a given matrix and particular operation type. The analysis meta data can be cleared by rocsparse_bsrsv_clear().rocsparse_bsrsv_analysis
can share its meta data with rocsparse_sbsrsm_analysis(), rocsparse_dbsrsm_analysis(), rocsparse_cbsrsm_analysis(), rocsparse_zbsrsm_analysis(), rocsparse_sbsrilu0_analysis(), rocsparse_dbsrilu0_analysis(), rocsparse_cbsrilu0_analysis(), rocsparse_zbsrilu0_analysis(), rocsparse_sbsric0_analysis(), rocsparse_dbsric0_analysis(), rocsparse_cbsric0_analysis() and rocsparse_zbsric0_analysis(). Selecting rocsparse_analysis_policy_reuse policy can greatly improve computation performance of meta data. However, the user need to make sure that the sparsity pattern remains unchanged. If this cannot be assured, rocsparse_analysis_policy_force has to be used.Note
If the matrix sparsity pattern changes, the gathered information will become invalid.
Note
This function is blocking with respect to the host.
Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
dir – [in] matrix storage of BSR blocks.
trans – [in] matrix operation type.
mb – [in] number of block rows of the sparse BSR matrix.
nnzb – [in] number of non-zero blocks of the sparse BSR matrix.
descr – [in] descriptor of the sparse BSR matrix.
bsr_val – [in] array of
nnzb
blocks of the sparse BSR matrix.bsr_row_ptr – [in] array of
mb+1
elements that point to the start of every block row of the sparse BSR matrix.bsr_col_ind – [in] array of
nnz
containing the block column indices of the sparse BSR matrix.block_dim – [in] block dimension of the sparse BSR matrix.
info – [out] structure that holds the information collected during the analysis step.
analysis – [in] rocsparse_analysis_policy_reuse or rocsparse_analysis_policy_force.
solve – [in] rocsparse_solve_policy_auto.
temp_buffer – [in] temporary storage buffer allocated by the user.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
mb
,nnzb
orblock_dim
is invalid.rocsparse_status_invalid_pointer –
descr
,bsr_row_ptr
,bsr_col_ind
,info
ortemp_buffer
pointer is invalid.rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_not_implemented –
trans
== rocsparse_operation_conjugate_transpose or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_bsrsv_solve()#
-
rocsparse_status rocsparse_sbsrsv_solve(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nnzb, const float *alpha, const rocsparse_mat_descr descr, const float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, const float *x, float *y, rocsparse_solve_policy policy, void *temp_buffer)#
-
rocsparse_status rocsparse_dbsrsv_solve(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nnzb, const double *alpha, const rocsparse_mat_descr descr, const double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, const double *x, double *y, rocsparse_solve_policy policy, void *temp_buffer)#
-
rocsparse_status rocsparse_cbsrsv_solve(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_float_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, const rocsparse_float_complex *x, rocsparse_float_complex *y, rocsparse_solve_policy policy, void *temp_buffer)#
-
rocsparse_status rocsparse_zbsrsv_solve(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_double_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, const rocsparse_double_complex *x, rocsparse_double_complex *y, rocsparse_solve_policy policy, void *temp_buffer)#
Sparse triangular solve using BSR storage format.
rocsparse_bsrsv_solve
solves a sparse triangular linear system of a sparse \(m \times m\) matrix, defined in BSR storage format, a dense solution vector \(y\) and the right-hand side \(x\) that is multiplied by \(\alpha\), such that\[ op(A) \cdot y = \alpha \cdot x, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == rocsparse_operation_none} \\ A^T, & \text{if trans == rocsparse_operation_transpose} \\ A^H, & \text{if trans == rocsparse_operation_conjugate_transpose} \end{array} \right. \end{split}\]rocsparse_bsrsv_solve
requires a user allocated temporary buffer. Its size is returned by rocsparse_sbsrsv_buffer_size(), rocsparse_dbsrsv_buffer_size(), rocsparse_cbsrsv_buffer_size() or rocsparse_zbsrsv_buffer_size(). Furthermore, analysis meta data is required. It can be obtained by rocsparse_sbsrsv_analysis(), rocsparse_dbsrsv_analysis(), rocsparse_cbsrsv_analysis() or rocsparse_zbsrsv_analysis().rocsparse_bsrsv_solve
reports the first zero pivot (either numerical or structural zero). The zero pivot status can be checked calling rocsparse_bsrsv_zero_pivot(). If rocsparse_diag_type == rocsparse_diag_type_unit, no zero pivot will be reported, even if \(A_{j,j} = 0\) for some \(j\).- Example
Consider the lower triangular \(m \times m\) matrix \(L\), stored in BSR storage format with unit diagonal. The following example solves \(L \cdot y = x\).
// Create rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); // Create matrix descriptor rocsparse_mat_descr descr; rocsparse_create_mat_descr(&descr); rocsparse_set_mat_fill_mode(descr, rocsparse_fill_mode_lower); rocsparse_set_mat_diag_type(descr, rocsparse_diag_type_unit); // Create matrix info structure rocsparse_mat_info info; rocsparse_create_mat_info(&info); // Obtain required buffer size size_t buffer_size; rocsparse_dbsrsv_buffer_size(handle, rocsparse_direction_column, rocsparse_operation_none, mb, nnzb, descr, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, &buffer_size); // Allocate temporary buffer void* temp_buffer; hipMalloc(&temp_buffer, buffer_size); // Perform analysis step rocsparse_dbsrsv_analysis(handle, rocsparse_direction_column, rocsparse_operation_none, mb, nnzb, descr, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); // Solve Ly = x rocsparse_dbsrsv_solve(handle, rocsparse_direction_column, rocsparse_operation_none, mb, nnzb, &alpha, descr, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, x, y, rocsparse_solve_policy_auto, temp_buffer); // No zero pivot should be found, with L having unit diagonal // Clean up hipFree(temp_buffer); rocsparse_destroy_mat_info(info); rocsparse_destroy_mat_descr(descr); rocsparse_destroy_handle(handle);
Note
The sparse BSR matrix has to be sorted.
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
Currently, only
trans
== rocsparse_operation_none andtrans
== rocsparse_operation_transpose is supported.Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
dir – [in] matrix storage of BSR blocks.
trans – [in] matrix operation type.
mb – [in] number of block rows of the sparse BSR matrix.
nnzb – [in] number of non-zero blocks of the sparse BSR matrix.
alpha – [in] scalar \(\alpha\).
descr – [in] descriptor of the sparse BSR matrix.
bsr_val – [in] array of
nnzb
blocks of the sparse BSR matrix.bsr_row_ptr – [in] array of
mb+1
elements that point to the start of every block row of the sparse BSR matrix.bsr_col_ind – [in] array of
nnz
containing the block column indices of the sparse BSR matrix.block_dim – [in] block dimension of the sparse BSR matrix.
info – [in] structure that holds the information collected during the analysis step.
x – [in] array of
m
elements, holding the right-hand side.y – [out] array of
m
elements, holding the solution.policy – [in] rocsparse_solve_policy_auto.
temp_buffer – [in] temporary storage buffer allocated by the user.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
mb
,nnzb
orblock_dim
is invalid.rocsparse_status_invalid_pointer –
descr
,alpha
,bsr_val
,bsr_row_ptr
,bsr_col_ind
,x
ory
pointer is invalid.rocsparse_status_arch_mismatch – the device is not supported.
rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_not_implemented –
trans
== rocsparse_operation_conjugate_transpose or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_bsrsv_clear()#
-
rocsparse_status rocsparse_bsrsv_clear(rocsparse_handle handle, rocsparse_mat_info info)#
Sparse triangular solve using BSR storage format.
rocsparse_bsrsv_clear
deallocates all memory that was allocated by rocsparse_sbsrsv_analysis(), rocsparse_dbsrsv_analysis(), rocsparse_cbsrsv_analysis() or rocsparse_zbsrsv_analysis(). This is especially useful, if memory is an issue and the analysis data is not required for further computation, e.g. when switching to another sparse matrix format. Callingrocsparse_bsrsv_clear
is optional. All allocated resources will be cleared, when the opaque rocsparse_mat_info struct is destroyed using rocsparse_destroy_mat_info().Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
info – [inout] structure that holds the information collected during the analysis step.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_pointer –
info
pointer is invalid.rocsparse_status_memory_error – the buffer holding the meta data could not be deallocated.
rocsparse_status_internal_error – an internal error occurred.
rocsparse_coomv()#
-
rocsparse_status rocsparse_scoomv(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, const float *alpha, const rocsparse_mat_descr descr, const float *coo_val, const rocsparse_int *coo_row_ind, const rocsparse_int *coo_col_ind, const float *x, const float *beta, float *y)#
-
rocsparse_status rocsparse_dcoomv(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, const double *alpha, const rocsparse_mat_descr descr, const double *coo_val, const rocsparse_int *coo_row_ind, const rocsparse_int *coo_col_ind, const double *x, const double *beta, double *y)#
-
rocsparse_status rocsparse_ccoomv(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, const rocsparse_float_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_float_complex *coo_val, const rocsparse_int *coo_row_ind, const rocsparse_int *coo_col_ind, const rocsparse_float_complex *x, const rocsparse_float_complex *beta, rocsparse_float_complex *y)#
-
rocsparse_status rocsparse_zcoomv(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, const rocsparse_double_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_double_complex *coo_val, const rocsparse_int *coo_row_ind, const rocsparse_int *coo_col_ind, const rocsparse_double_complex *x, const rocsparse_double_complex *beta, rocsparse_double_complex *y)#
Sparse matrix vector multiplication using COO storage format.
rocsparse_coomv
multiplies the scalar \(\alpha\) with a sparse \(m \times n\) matrix, defined in COO storage format, and the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == rocsparse_operation_none} \\ A^T, & \text{if trans == rocsparse_operation_transpose} \\ A^H, & \text{if trans == rocsparse_operation_conjugate_transpose} \end{array} \right. \end{split}\]The COO matrix has to be sorted by row indices. This can be achieved by using rocsparse_coosort_by_row().
for(i = 0; i < m; ++i) { y[i] = beta * y[i]; } for(i = 0; i < nnz; ++i) { y[coo_row_ind[i]] += alpha * coo_val[i] * x[coo_col_ind[i]]; }
- Example
This example performs a sparse matrix vector multiplication in COO format.
// rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); // A sparse matrix // 1 0 3 4 // 0 0 5 1 // 0 2 0 0 // 4 0 0 8 rocsparse_int hArow[8] = {0, 0, 0, 1, 1, 2, 3, 3}; rocsparse_int hAcol[8] = {0, 2, 3, 2, 3, 1, 0, 3}; double hAval[8] = {1.0, 3.0, 4.0, 5.0, 1.0, 2.0, 4.0, 8.0}; rocsparse_int m = 4; rocsparse_int n = 4; rocsparse_int nnz = 8; double halpha = 1.0; double hbeta = 0.0; double hx[4] = {1.0, 2.0, 3.0, 4.0}; // Matrix descriptor rocsparse_mat_descr descrA; rocsparse_create_mat_descr(&descrA); // Offload data to device rocsparse_int* dArow = NULL; rocsparse_int* dAcol = NULL; double* dAval = NULL; double* dx = NULL; double* dy = NULL; hipMalloc((void**)&dArow, sizeof(rocsparse_int) * nnz); hipMalloc((void**)&dAcol, sizeof(rocsparse_int) * nnz); hipMalloc((void**)&dAval, sizeof(double) * nnz); hipMalloc((void**)&dx, sizeof(double) * n); hipMalloc((void**)&dy, sizeof(double) * m); hipMemcpy(dArow, hArow, sizeof(rocsparse_int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dAcol, hAcol, sizeof(rocsparse_int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dAval, hAval, sizeof(double) * nnz, hipMemcpyHostToDevice); hipMemcpy(dx, hx, sizeof(double) * n, hipMemcpyHostToDevice); // Call rocsparse coomv rocsparse_dcoomv(handle, rocsparse_operation_none, m, n, nnz, &halpha, descrA, dAval, dArow, dAcol, dx, &hbeta, dy); // Copy back to host double hy[4]; hipMemcpy(hy, dy, sizeof(double) * m, hipMemcpyDeviceToHost); // Clear up on device hipFree(dArow); hipFree(dAcol); hipFree(dAval); hipFree(dx); hipFree(dy); rocsparse_destroy_mat_descr(descrA); rocsparse_destroy_handle(handle);
Note
This function does not produce deterministic results when A is transposed.
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
trans – [in] matrix operation type.
m – [in] number of rows of the sparse COO matrix.
n – [in] number of columns of the sparse COO matrix.
nnz – [in] number of non-zero entries of the sparse COO matrix.
alpha – [in] scalar \(\alpha\).
descr – [in] descriptor of the sparse COO matrix. Currently, only rocsparse_matrix_type_general is supported.
coo_val – [in] array of
nnz
elements of the sparse COO matrix.coo_row_ind – [in] array of
nnz
elements containing the row indices of the sparse COO matrix.coo_col_ind – [in] array of
nnz
elements containing the column indices of the sparse COO matrix.x – [in] array of
n
elements ( \(op(A) = A\)) orm
elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).beta – [in] scalar \(\beta\).
y – [inout] array of
m
elements ( \(op(A) = A\)) orn
elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
m
,n
ornnz
is invalid.rocsparse_status_invalid_pointer –
descr
,alpha
,coo_val
,coo_row_ind
,coo_col_ind
,x
,beta
ory
pointer is invalid.rocsparse_status_arch_mismatch – the device is not supported.
rocsparse_status_not_implemented – rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_csrmv_analysis()#
-
rocsparse_status rocsparse_scsrmv_analysis(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, const rocsparse_mat_descr descr, const float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info)#
-
rocsparse_status rocsparse_dcsrmv_analysis(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, const rocsparse_mat_descr descr, const double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info)#
-
rocsparse_status rocsparse_ccsrmv_analysis(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info)#
-
rocsparse_status rocsparse_zcsrmv_analysis(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info)#
Sparse matrix vector multiplication using CSR storage format.
rocsparse_csrmv_analysis
performs the analysis step for rocsparse_scsrmv(), rocsparse_dcsrmv(), rocsparse_ccsrmv() and rocsparse_zcsrmv(). It is expected that this function will be executed only once for a given matrix and particular operation type. The gathered analysis meta data can be cleared by rocsparse_csrmv_clear().Note
If the matrix sparsity pattern changes, the gathered information will become invalid.
Note
This function is blocking with respect to the host.
Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
trans – [in] matrix operation type.
m – [in] number of rows of the sparse CSR matrix.
n – [in] number of columns of the sparse CSR matrix.
nnz – [in] number of non-zero entries of the sparse CSR matrix.
descr – [in] descriptor of the sparse CSR matrix.
csr_val – [in] array of
nnz
elements of the sparse CSR matrix.csr_row_ptr – [in] array of
m+1
elements that point to the start of every row of the sparse CSR matrix.csr_col_ind – [in] array of
nnz
elements containing the column indices of the sparse CSR matrix.info – [out] structure that holds the information collected during the analysis step.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
m
,n
ornnz
is invalid.rocsparse_status_invalid_pointer –
descr
,csr_val
,csr_row_ptr
,csr_col_ind
orinfo
pointer is invalid.rocsparse_status_memory_error – the buffer for the gathered information could not be allocated.
rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_not_implemented – if rocsparse_matrix_type is not one of rocsparse_matrix_type_general, rocsparse_matrix_type_symmetric, or rocsparse_matrix_type_triangular.
rocsparse_csrmv()#
-
rocsparse_status rocsparse_scsrmv(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, const float *alpha, const rocsparse_mat_descr descr, const float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, const float *x, const float *beta, float *y)#
-
rocsparse_status rocsparse_dcsrmv(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, const double *alpha, const rocsparse_mat_descr descr, const double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, const double *x, const double *beta, double *y)#
-
rocsparse_status rocsparse_ccsrmv(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, const rocsparse_float_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, const rocsparse_float_complex *x, const rocsparse_float_complex *beta, rocsparse_float_complex *y)#
-
rocsparse_status rocsparse_zcsrmv(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, const rocsparse_double_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, const rocsparse_double_complex *x, const rocsparse_double_complex *beta, rocsparse_double_complex *y)#
Sparse matrix vector multiplication using CSR storage format.
rocsparse_csrmv
multiplies the scalar \(\alpha\) with a sparse \(m \times n\) matrix, defined in CSR storage format, and the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == rocsparse_operation_none} \\ A^T, & \text{if trans == rocsparse_operation_transpose} \\ A^H, & \text{if trans == rocsparse_operation_conjugate_transpose} \end{array} \right. \end{split}\]The
info
parameter is optional and contains information collected by rocsparse_scsrmv_analysis(), rocsparse_dcsrmv_analysis(), rocsparse_ccsrmv_analysis() or rocsparse_zcsrmv_analysis(). If present, the information will be used to speed up thecsrmv
computation. Ifinfo
==NULL
, generalcsrmv
routine will be used instead.for(i = 0; i < m; ++i) { y[i] = beta * y[i]; for(j = csr_row_ptr[i]; j < csr_row_ptr[i + 1]; ++j) { y[i] = y[i] + alpha * csr_val[j] * x[csr_col_ind[j]]; } }
- Example
This example performs a sparse matrix vector multiplication in CSR format using additional meta data to improve performance.
// Create matrix info structure rocsparse_mat_info info; rocsparse_create_mat_info(&info); // Perform analysis step to obtain meta data rocsparse_scsrmv_analysis(handle, rocsparse_operation_none, m, n, nnz, descr, csr_val, csr_row_ptr, csr_col_ind, info); // Compute y = Ax rocsparse_scsrmv(handle, rocsparse_operation_none, m, n, nnz, &alpha, descr, csr_val, csr_row_ptr, csr_col_ind, info, x, &beta, y); // Do more work // ... // Clean up rocsparse_destroy_mat_info(info);
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
trans – [in] matrix operation type.
m – [in] number of rows of the sparse CSR matrix.
n – [in] number of columns of the sparse CSR matrix.
nnz – [in] number of non-zero entries of the sparse CSR matrix.
alpha – [in] scalar \(\alpha\).
descr – [in] descriptor of the sparse CSR matrix. Currently, only rocsparse_matrix_type_general is supported.
csr_val – [in] array of
nnz
elements of the sparse CSR matrix.csr_row_ptr – [in] array of
m+1
elements that point to the start of every row of the sparse CSR matrix.csr_col_ind – [in] array of
nnz
elements containing the column indices of the sparse CSR matrix.info – [in] information collected by rocsparse_scsrmv_analysis(), rocsparse_dcsrmv_analysis(), rocsparse_ccsrmv_analysis() or rocsparse_dcsrmv_analysis(), can be
NULL
if no information is available.x – [in] array of
n
elements ( \(op(A) == A\)) orm
elements ( \(op(A) == A^T\) or \(op(A) == A^H\)).beta – [in] scalar \(\beta\).
y – [inout] array of
m
elements ( \(op(A) == A\)) orn
elements ( \(op(A) == A^T\) or \(op(A) == A^H\)).
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
m
,n
ornnz
is invalid.rocsparse_status_invalid_pointer –
descr
,alpha
,csr_val
,csr_row_ptr
,csr_col_ind
,x
,beta
ory
pointer is invalid.rocsparse_status_arch_mismatch – the device is not supported.
rocsparse_status_not_implemented – rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_csrmv_analysis_clear()#
-
rocsparse_status rocsparse_csrmv_clear(rocsparse_handle handle, rocsparse_mat_info info)#
Sparse matrix vector multiplication using CSR storage format.
rocsparse_csrmv_clear
deallocates all memory that was allocated by rocsparse_scsrmv_analysis(), rocsparse_dcsrmv_analysis(), rocsparse_ccsrmv_analysis() or rocsparse_zcsrmv_analysis(). This is especially useful, if memory is an issue and the analysis data is not required anymore for further computation, e.g. when switching to another sparse matrix format.Note
Calling
rocsparse_csrmv_clear
is optional. All allocated resources will be cleared, when the opaque rocsparse_mat_info struct is destroyed using rocsparse_destroy_mat_info().Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
info – [inout] structure that holds the information collected during analysis step.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_pointer –
info
pointer is invalid.rocsparse_status_memory_error – the buffer for the gathered information could not be deallocated.
rocsparse_status_internal_error – an internal error occurred.
rocsparse_csrsv_zero_pivot()#
-
rocsparse_status rocsparse_csrsv_zero_pivot(rocsparse_handle handle, const rocsparse_mat_descr descr, rocsparse_mat_info info, rocsparse_int *position)#
Sparse triangular solve using CSR storage format.
rocsparse_csrsv_zero_pivot
returns rocsparse_status_zero_pivot, if either a structural or numerical zero has been found during rocsparse_scsrsv_solve(), rocsparse_dcsrsv_solve(), rocsparse_ccsrsv_solve() or rocsparse_zcsrsv_solve() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored inposition
, using same index base as the CSR matrix.position
can be in host or device memory. If no zero pivot has been found,position
is set to -1 and rocsparse_status_success is returned instead.Note
rocsparse_csrsv_zero_pivot
is a blocking function. It might influence performance negatively.Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
descr – [in] descriptor of the sparse CSR matrix.
info – [in] structure that holds the information collected during the analysis step.
position – [inout] pointer to zero pivot \(j\), can be in host or device memory.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_pointer –
info
orposition
pointer is invalid.rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_zero_pivot – zero pivot has been found.
rocsparse_csrsv_buffer_size()#
-
rocsparse_status rocsparse_scsrsv_buffer_size(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
-
rocsparse_status rocsparse_dcsrsv_buffer_size(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
-
rocsparse_status rocsparse_ccsrsv_buffer_size(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
-
rocsparse_status rocsparse_zcsrsv_buffer_size(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
Sparse triangular solve using CSR storage format.
rocsparse_csrsv_buffer_size
returns the size of the temporary storage buffer that is required by rocsparse_scsrsv_analysis(), rocsparse_dcsrsv_analysis(), rocsparse_ccsrsv_analysis(), rocsparse_zcsrsv_analysis(), rocsparse_scsrsv_solve(), rocsparse_dcsrsv_solve(), rocsparse_ccsrsv_solve() and rocsparse_zcsrsv_solve(). The temporary storage buffer must be allocated by the user. The size of the temporary storage buffer is identical to the size returned by rocsparse_scsrilu0_buffer_size(), rocsparse_dcsrilu0_buffer_size(), rocsparse_ccsrilu0_buffer_size() and rocsparse_zcsrilu0_buffer_size() if the matrix sparsity pattern is identical. The user allocated buffer can thus be shared between subsequent calls to those functions.Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
trans – [in] matrix operation type.
m – [in] number of rows of the sparse CSR matrix.
nnz – [in] number of non-zero entries of the sparse CSR matrix.
descr – [in] descriptor of the sparse CSR matrix.
csr_val – [in] array of
nnz
elements of the sparse CSR matrix.csr_row_ptr – [in] array of
m+1
elements that point to the start of every row of the sparse CSR matrix.csr_col_ind – [in] array of
nnz
elements containing the column indices of the sparse CSR matrix.info – [out] structure that holds the information collected during the analysis step.
buffer_size – [out] number of bytes of the temporary storage buffer required by rocsparse_scsrsv_analysis(), rocsparse_dcsrsv_analysis(), rocsparse_ccsrsv_analysis(), rocsparse_zcsrsv_analysis(), rocsparse_scsrsv_solve(), rocsparse_dcsrsv_solve(), rocsparse_ccsrsv_solve() and rocsparse_zcsrsv_solve().
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
m
ornnz
is invalid.rocsparse_status_invalid_pointer –
descr
,csr_val
,csr_row_ptr
,csr_col_ind
,info
orbuffer_size
pointer is invalid.rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_not_implemented –
trans
== rocsparse_operation_conjugate_transpose or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_csrsv_analysis()#
-
rocsparse_status rocsparse_scsrsv_analysis(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
-
rocsparse_status rocsparse_dcsrsv_analysis(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
-
rocsparse_status rocsparse_ccsrsv_analysis(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
-
rocsparse_status rocsparse_zcsrsv_analysis(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
Sparse triangular solve using CSR storage format.
rocsparse_csrsv_analysis
performs the analysis step for rocsparse_scsrsv_solve(), rocsparse_dcsrsv_solve(), rocsparse_ccsrsv_solve() and rocsparse_zcsrsv_solve(). It is expected that this function will be executed only once for a given matrix and particular operation type. The analysis meta data can be cleared by rocsparse_csrsv_clear().rocsparse_csrsv_analysis
can share its meta data with rocsparse_scsrsm_analysis(), rocsparse_dcsrsm_analysis(), rocsparse_ccsrsm_analysis(), rocsparse_zcsrsm_analysis(), rocsparse_scsrilu0_analysis(), rocsparse_dcsrilu0_analysis(), rocsparse_ccsrilu0_analysis(), rocsparse_zcsrilu0_analysis(), rocsparse_scsric0_analysis(), rocsparse_dcsric0_analysis(), rocsparse_ccsric0_analysis() and rocsparse_zcsric0_analysis(). Selecting rocsparse_analysis_policy_reuse policy can greatly improve computation performance of meta data. However, the user need to make sure that the sparsity pattern remains unchanged. If this cannot be assured, rocsparse_analysis_policy_force has to be used.Note
If the matrix sparsity pattern changes, the gathered information will become invalid.
Note
This function is blocking with respect to the host.
Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
trans – [in] matrix operation type.
m – [in] number of rows of the sparse CSR matrix.
nnz – [in] number of non-zero entries of the sparse CSR matrix.
descr – [in] descriptor of the sparse CSR matrix.
csr_val – [in] array of
nnz
elements of the sparse CSR matrix.csr_row_ptr – [in] array of
m+1
elements that point to the start of every row of the sparse CSR matrix.csr_col_ind – [in] array of
nnz
elements containing the column indices of the sparse CSR matrix.info – [out] structure that holds the information collected during the analysis step.
analysis – [in] rocsparse_analysis_policy_reuse or rocsparse_analysis_policy_force.
solve – [in] rocsparse_solve_policy_auto.
temp_buffer – [in] temporary storage buffer allocated by the user.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
m
ornnz
is invalid.rocsparse_status_invalid_pointer –
descr
,csr_row_ptr
,csr_col_ind
,info
ortemp_buffer
pointer is invalid.rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_not_implemented –
trans
== rocsparse_operation_conjugate_transpose or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_csrsv_solve()#
-
rocsparse_status rocsparse_scsrsv_solve(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const float *alpha, const rocsparse_mat_descr descr, const float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, const float *x, float *y, rocsparse_solve_policy policy, void *temp_buffer)#
-
rocsparse_status rocsparse_dcsrsv_solve(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const double *alpha, const rocsparse_mat_descr descr, const double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, const double *x, double *y, rocsparse_solve_policy policy, void *temp_buffer)#
-
rocsparse_status rocsparse_ccsrsv_solve(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_float_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, const rocsparse_float_complex *x, rocsparse_float_complex *y, rocsparse_solve_policy policy, void *temp_buffer)#
-
rocsparse_status rocsparse_zcsrsv_solve(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_double_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, const rocsparse_double_complex *x, rocsparse_double_complex *y, rocsparse_solve_policy policy, void *temp_buffer)#
Sparse triangular solve using CSR storage format.
rocsparse_csrsv_solve
solves a sparse triangular linear system of a sparse \(m \times m\) matrix, defined in CSR storage format, a dense solution vector \(y\) and the right-hand side \(x\) that is multiplied by \(\alpha\), such that\[ op(A) \cdot y = \alpha \cdot x, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == rocsparse_operation_none} \\ A^T, & \text{if trans == rocsparse_operation_transpose} \\ A^H, & \text{if trans == rocsparse_operation_conjugate_transpose} \end{array} \right. \end{split}\]rocsparse_csrsv_solve
requires a user allocated temporary buffer. Its size is returned by rocsparse_scsrsv_buffer_size(), rocsparse_dcsrsv_buffer_size(), rocsparse_ccsrsv_buffer_size() or rocsparse_zcsrsv_buffer_size(). Furthermore, analysis meta data is required. It can be obtained by rocsparse_scsrsv_analysis(), rocsparse_dcsrsv_analysis(), rocsparse_ccsrsv_analysis() or rocsparse_zcsrsv_analysis().rocsparse_csrsv_solve
reports the first zero pivot (either numerical or structural zero). The zero pivot status can be checked calling rocsparse_csrsv_zero_pivot(). If rocsparse_diag_type == rocsparse_diag_type_unit, no zero pivot will be reported, even if \(A_{j,j} = 0\) for some \(j\).- Example
Consider the lower triangular \(m \times m\) matrix \(L\), stored in CSR storage format with unit diagonal. The following example solves \(L \cdot y = x\).
// Create rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); // Create matrix descriptor rocsparse_mat_descr descr; rocsparse_create_mat_descr(&descr); rocsparse_set_mat_fill_mode(descr, rocsparse_fill_mode_lower); rocsparse_set_mat_diag_type(descr, rocsparse_diag_type_unit); // Create matrix info structure rocsparse_mat_info info; rocsparse_create_mat_info(&info); // Obtain required buffer size size_t buffer_size; rocsparse_dcsrsv_buffer_size(handle, rocsparse_operation_none, m, nnz, descr, csr_val, csr_row_ptr, csr_col_ind, info, &buffer_size); // Allocate temporary buffer void* temp_buffer; hipMalloc(&temp_buffer, buffer_size); // Perform analysis step rocsparse_dcsrsv_analysis(handle, rocsparse_operation_none, m, nnz, descr, csr_val, csr_row_ptr, csr_col_ind, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); // Solve Ly = x rocsparse_dcsrsv_solve(handle, rocsparse_operation_none, m, nnz, &alpha, descr, csr_val, csr_row_ptr, csr_col_ind, info, x, y, rocsparse_solve_policy_auto, temp_buffer); // No zero pivot should be found, with L having unit diagonal // Clean up hipFree(temp_buffer); rocsparse_destroy_mat_info(info); rocsparse_destroy_mat_descr(descr); rocsparse_destroy_handle(handle);
Note
The sparse CSR matrix has to be sorted. This can be achieved by calling rocsparse_csrsort().
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
Currently, only
trans
== rocsparse_operation_none andtrans
== rocsparse_operation_transpose is supported.Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
trans – [in] matrix operation type.
m – [in] number of rows of the sparse CSR matrix.
nnz – [in] number of non-zero entries of the sparse CSR matrix.
alpha – [in] scalar \(\alpha\).
descr – [in] descriptor of the sparse CSR matrix.
csr_val – [in] array of
nnz
elements of the sparse CSR matrix.csr_row_ptr – [in] array of
m+1
elements that point to the start of every row of the sparse CSR matrix.csr_col_ind – [in] array of
nnz
elements containing the column indices of the sparse CSR matrix.info – [in] structure that holds the information collected during the analysis step.
x – [in] array of
m
elements, holding the right-hand side.y – [out] array of
m
elements, holding the solution.policy – [in] rocsparse_solve_policy_auto.
temp_buffer – [in] temporary storage buffer allocated by the user.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
m
ornnz
is invalid.rocsparse_status_invalid_pointer –
descr
,alpha
,csr_val
,csr_row_ptr
,csr_col_ind
,x
ory
pointer is invalid.rocsparse_status_arch_mismatch – the device is not supported.
rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_not_implemented –
trans
== rocsparse_operation_conjugate_transpose or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_csrsv_clear()#
-
rocsparse_status rocsparse_csrsv_clear(rocsparse_handle handle, const rocsparse_mat_descr descr, rocsparse_mat_info info)#
Sparse triangular solve using CSR storage format.
rocsparse_csrsv_clear
deallocates all memory that was allocated by rocsparse_scsrsv_analysis(), rocsparse_dcsrsv_analysis(), rocsparse_ccsrsv_analysis() or rocsparse_zcsrsv_analysis(). This is especially useful, if memory is an issue and the analysis data is not required for further computation, e.g. when switching to another sparse matrix format. Callingrocsparse_csrsv_clear
is optional. All allocated resources will be cleared, when the opaque rocsparse_mat_info struct is destroyed using rocsparse_destroy_mat_info().Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
descr – [in] descriptor of the sparse CSR matrix.
info – [inout] structure that holds the information collected during the analysis step.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_pointer –
info
pointer is invalid.rocsparse_status_memory_error – the buffer holding the meta data could not be deallocated.
rocsparse_status_internal_error – an internal error occurred.
rocsparse_csritsv_zero_pivot()#
-
rocsparse_status rocsparse_csritsv_zero_pivot(rocsparse_handle handle, const rocsparse_mat_descr descr, rocsparse_mat_info info, rocsparse_int *position)#
Sparse iterative triangular solve using CSR storage format.
rocsparse_csritsv_zero_pivot
returns rocsparse_status_zero_pivot, if either a structural or numerical zero has been found during rocsparse_csritsv_solve() and or rocsparse_csritsv_analysis(), execution. The first zero pivot \(j\) at \(A_{j,j}\) is stored inposition
, using same index base as the CSR matrix.position
can be in host or device memory. If no zero pivot has been found,position
is set to -1 and rocsparse_status_success is returned instead.Note
rocsparse_csritsv_zero_pivot
is a blocking function. It might influence performance negatively.Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
descr – [in] descriptor of the sparse CSR matrix.
info – [in] structure that holds the information collected during the analysis step.
position – [inout] pointer to zero pivot \(j\), can be in host or device memory.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_pointer –
info
orposition
pointer is invalid.rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_zero_pivot – zero pivot has been found.
rocsparse_csritsv_buffer_size()#
-
rocsparse_status rocsparse_scsritsv_buffer_size(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
-
rocsparse_status rocsparse_dcsritsv_buffer_size(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
-
rocsparse_status rocsparse_ccsritsv_buffer_size(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
-
rocsparse_status rocsparse_zcsritsv_buffer_size(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
Sparse iterative triangular solve using CSR storage format.
rocsparse_csritsv_buffer_size
returns the size of the temporary storage buffer that is required by rocsparse_scsritsv_analysis(), rocsparse_dcsritsv_analysis(), rocsparse_ccsritsv_analysis(), rocsparse_zcsritsv_analysis(), rocsparse_scsritsv_solve(), rocsparse_dcsritsv_solve(), rocsparse_ccsritsv_solve() and rocsparse_zcsritsv_solve(). The temporary storage buffer must be allocated by the user.Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
trans – [in] matrix operation type.
m – [in] number of rows of the sparse CSR matrix.
nnz – [in] number of non-zero entries of the sparse CSR matrix.
descr – [in] descriptor of the sparse CSR matrix.
csr_val – [in] array of
nnz
elements of the sparse CSR matrix.csr_row_ptr – [in] array of
m+1
elements that point to the start of every row of the sparse CSR matrix.csr_col_ind – [in] array of
nnz
elements containing the column indices of the sparse CSR matrix.info – [out] structure that holds the information collected during the analysis step.
buffer_size – [out] number of bytes of the temporary storage buffer required by rocsparse_scsritsv_analysis(), rocsparse_dcsritsv_analysis(), rocsparse_ccsritsv_analysis(), rocsparse_zcsritsv_analysis(), rocsparse_scsritsv_solve(), rocsparse_dcsritsv_solve(), rocsparse_ccsritsv_solve() and rocsparse_zcsritsv_solve().
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
m
ornnz
is invalid.rocsparse_status_invalid_pointer –
descr
,csr_val
,csr_row_ptr
,csr_col_ind
,info
orbuffer_size
pointer is invalid.rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_not_implemented – rocsparse_matrix_type != rocsparse_matrix_type_general and rocsparse_matrix_type != rocsparse_matrix_type_triangular.
rocsparse_csritsv_analysis()#
-
rocsparse_status rocsparse_scsritsv_analysis(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
-
rocsparse_status rocsparse_dcsritsv_analysis(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
-
rocsparse_status rocsparse_ccsritsv_analysis(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
-
rocsparse_status rocsparse_zcsritsv_analysis(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
Sparse iterative triangular solve using CSR storage format.
rocsparse_csritsv_analysis
performs the analysis step for rocsparse_scsritsv_solve(), rocsparse_dcsritsv_solve(), rocsparse_ccsritsv_solve() and rocsparse_zcsritsv_solve(). It is expected that this function will be executed only once for a given matrix and particular operation type. The analysis meta data can be cleared by rocsparse_csritsv_clear().Selecting rocsparse_analysis_policy_reuse policy can greatly improve computation performance of meta data. However, the user need to make sure that the sparsity pattern remains unchanged. If this cannot be assured, rocsparse_analysis_policy_force has to be used.
Note
If the matrix sparsity pattern changes, the gathered information will become invalid.
Note
This function is blocking with respect to the host.
Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
trans – [in] matrix operation type.
m – [in] number of rows of the sparse CSR matrix.
nnz – [in] number of non-zero entries of the sparse CSR matrix.
descr – [in] descriptor of the sparse CSR matrix.
csr_val – [in] array of
nnz
elements of the sparse CSR matrix.csr_row_ptr – [in] array of
m+1
elements that point to the start of every row of the sparse CSR matrix.csr_col_ind – [in] array of
nnz
elements containing the column indices of the sparse CSR matrix.info – [out] structure that holds the information collected during the analysis step.
analysis – [in] rocsparse_analysis_policy_reuse or rocsparse_analysis_policy_force.
solve – [in] rocsparse_solve_policy_auto.
temp_buffer – [in] temporary storage buffer allocated by the user.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
m
ornnz
is invalid.rocsparse_status_invalid_pointer –
descr
,csr_row_ptr
,csr_col_ind
,info
ortemp_buffer
pointer is invalid.rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_not_implemented – rocsparse_matrix_type != rocsparse_matrix_type_general and rocsparse_matrix_type != rocsparse_matrix_type_triangular.
rocsparse_csritsv_solve()#
-
rocsparse_status rocsparse_scsritsv_solve(rocsparse_handle handle, rocsparse_int *host_nmaxiter, const float *host_tol, float *host_history, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const float *alpha, const rocsparse_mat_descr descr, const float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, const float *x, float *y, rocsparse_solve_policy policy, void *temp_buffer)#
-
rocsparse_status rocsparse_dcsritsv_solve(rocsparse_handle handle, rocsparse_int *host_nmaxiter, const double *host_tol, double *host_history, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const double *alpha, const rocsparse_mat_descr descr, const double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, const double *x, double *y, rocsparse_solve_policy policy, void *temp_buffer)#
-
rocsparse_status rocsparse_ccsritsv_solve(rocsparse_handle handle, rocsparse_int *host_nmaxiter, const float *host_tol, float *host_history, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_float_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, const rocsparse_float_complex *x, rocsparse_float_complex *y, rocsparse_solve_policy policy, void *temp_buffer)#
-
rocsparse_status rocsparse_zcsritsv_solve(rocsparse_handle handle, rocsparse_int *host_nmaxiter, const double *host_tol, double *host_history, rocsparse_operation trans, rocsparse_int m, rocsparse_int nnz, const rocsparse_double_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, const rocsparse_double_complex *x, rocsparse_double_complex *y, rocsparse_solve_policy policy, void *temp_buffer)#
Sparse iterative triangular solve using CSR storage format.
rocsparse_csritsv_solve
solves iteratively with the use of the Jacobi method a sparse triangular linear system of a sparse \(m \times m\) matrix, defined in CSR storage format, a dense solution vector \(y\) and the right-hand side \(x\) that is multiplied by \(\alpha\), such that\[ op(A) \cdot y = \alpha \cdot x, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == rocsparse_operation_none} \\ A^T, & \text{if trans == rocsparse_operation_transpose} \\ A^H, & \text{if trans == rocsparse_operation_conjugate_transpose} \end{array} \right. \end{split}\]rocsparse_csritsv_solve
requires a user allocated temporary buffer. Its size is returned by rocsparse_scsritsv_buffer_size(), rocsparse_dcsritsv_buffer_size(), rocsparse_ccsritsv_buffer_size() or rocsparse_zcsritsv_buffer_size(). Furthermore, analysis meta data is required. It can be obtained by rocsparse_scsritsv_analysis(), rocsparse_dcsritsv_analysis(), rocsparse_ccsritsv_analysis() or rocsparse_zcsritsv_analysis().rocsparse_csritsv_solve
reports the first zero pivot (either numerical or structural zero). The zero pivot status can be checked calling rocsparse_csritsv_zero_pivot(). If rocsparse_diag_type == rocsparse_diag_type_unit, no zero pivot will be reported, even if \(A_{j,j} = 0\) for some \(j\).- Example
Consider the lower triangular \(m \times m\) matrix \(L\), stored in CSR storage format with unit diagonal. The following example solves \(L \cdot y = x\).
// Create rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); // Create matrix descriptor rocsparse_mat_descr descr; rocsparse_create_mat_descr(&descr); rocsparse_set_mat_fill_mode(descr, rocsparse_fill_mode_lower); rocsparse_set_mat_diag_type(descr, rocsparse_diag_type_unit); // Create matrix info structure rocsparse_mat_info info; rocsparse_create_mat_info(&info); // Obtain required buffer size size_t buffer_size; rocsparse_dcsritsv_buffer_size(handle, rocsparse_operation_none, m, nnz, descr, csr_val, csr_row_ptr, csr_col_ind, info, &buffer_size); // Allocate temporary buffer void* temp_buffer; hipMalloc(&temp_buffer, buffer_size); // Perform analysis step rocsparse_dcsritsv_analysis(handle, rocsparse_operation_none, m, nnz, descr, csr_val, csr_row_ptr, csr_col_ind, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); // Solve Ly = x rocsparse_int nmaxiter = 200; rocsparse_int maxiter = nmaxiter; tol = 1.0e-4; history[200]; rocsparse_dcsritsv_solve(handle, &maxiter, &tol, history, rocsparse_operation_none, m, nnz, &alpha, descr, csr_val, csr_row_ptr, csr_col_ind, info, x, y, rocsparse_solve_policy_auto, temp_buffer); if (maxiter < nmaxiter) {} // convergence else {} // non converged for (int i=0;i<maxiter;++i) printf("iter = %d, max residual=%e\n", iter, history[i]); // No zero pivot should be found, with L having unit diagonal // Clean up hipFree(temp_buffer); rocsparse_destroy_mat_info(info); rocsparse_destroy_mat_descr(descr); rocsparse_destroy_handle(handle);
Note
The sparse CSR matrix has to be sorted. This can be achieved by calling rocsparse_csrsort().
Note
This function is blocking with respect to the host.
Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
host_nmaxiter – [inout] maximum number of iteration on input and maximum number of iteration on output.
host_tol – [in] if the pointer is null then loop will execute
nmaxiter
[0] iterations.host_history – [out] (optional, record history)
trans – [in] matrix operation type.
m – [in] number of rows of the sparse CSR matrix.
nnz – [in] number of non-zero entries of the sparse CSR matrix.
alpha – [in] scalar \(\alpha\).
descr – [in] descriptor of the sparse CSR matrix.
csr_val – [in] array of
nnz
elements of the sparse CSR matrix.csr_row_ptr – [in] array of
m+1
elements that point to the start of every row of the sparse CSR matrix.csr_col_ind – [in] array of
nnz
elements containing the column indices of the sparse CSR matrix.info – [in] structure that holds the information collected during the analysis step.
x – [in] array of
m
elements, holding the right-hand side.y – [out] array of
m
elements, holding the solution.policy – [in] rocsparse_solve_policy_auto.
temp_buffer – [in] temporary storage buffer allocated by the user.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
m
ornnz
is invalid.rocsparse_status_invalid_pointer –
descr
,alpha
,csr_val
,csr_row_ptr
,csr_col_ind
,x
ory
pointer is invalid.rocsparse_status_arch_mismatch – the device is not supported.
rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_not_implemented – rocsparse_matrix_type != rocsparse_matrix_type_general and rocsparse_matrix_type != rocsparse_matrix_type_triangular.
rocsparse_csritsv_clear()#
-
rocsparse_status rocsparse_csritsv_clear(rocsparse_handle handle, const rocsparse_mat_descr descr, rocsparse_mat_info info)#
Sparse triangular solve using CSR storage format.
rocsparse_csritsv_clear
deallocates all memory that was allocated by rocsparse_scsritsv_analysis(), rocsparse_dcsritsv_analysis(), rocsparse_ccsritsv_analysis() or rocsparse_zcsritsv_analysis(). This is especially useful, if memory is an issue and the analysis data is not required for further computation, e.g. when switching to another sparse matrix format. Callingrocsparse_csritsv_clear
is optional. All allocated resources will be cleared, when the opaque rocsparse_mat_info struct is destroyed using rocsparse_destroy_mat_info().Note
This routine does not support execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
descr – [in] descriptor of the sparse CSR matrix.
info – [inout] structure that holds the information collected during the analysis step.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_pointer –
info
pointer is invalid.rocsparse_status_memory_error – the buffer holding the meta data could not be deallocated.
rocsparse_status_internal_error – an internal error occurred.
rocsparse_ellmv()#
-
rocsparse_status rocsparse_sellmv(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, const float *alpha, const rocsparse_mat_descr descr, const float *ell_val, const rocsparse_int *ell_col_ind, rocsparse_int ell_width, const float *x, const float *beta, float *y)#
-
rocsparse_status rocsparse_dellmv(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, const double *alpha, const rocsparse_mat_descr descr, const double *ell_val, const rocsparse_int *ell_col_ind, rocsparse_int ell_width, const double *x, const double *beta, double *y)#
-
rocsparse_status rocsparse_cellmv(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, const rocsparse_float_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_float_complex *ell_val, const rocsparse_int *ell_col_ind, rocsparse_int ell_width, const rocsparse_float_complex *x, const rocsparse_float_complex *beta, rocsparse_float_complex *y)#
-
rocsparse_status rocsparse_zellmv(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, const rocsparse_double_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_double_complex *ell_val, const rocsparse_int *ell_col_ind, rocsparse_int ell_width, const rocsparse_double_complex *x, const rocsparse_double_complex *beta, rocsparse_double_complex *y)#
Sparse matrix vector multiplication using ELL storage format.
rocsparse_ellmv
multiplies the scalar \(\alpha\) with a sparse \(m \times n\) matrix, defined in ELL storage format, and the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == rocsparse_operation_none} \\ A^T, & \text{if trans == rocsparse_operation_transpose} \\ A^H, & \text{if trans == rocsparse_operation_conjugate_transpose} \end{array} \right. \end{split}\]for(i = 0; i < m; ++i) { y[i] = beta * y[i]; for(p = 0; p < ell_width; ++p) { idx = p * m + i; if((ell_col_ind[idx] >= 0) && (ell_col_ind[idx] < n)) { y[i] = y[i] + alpha * ell_val[idx] * x[ell_col_ind[idx]]; } } }
- Example
This example performs a sparse matrix vector multiplication in ELL format. It also shows how to convert from CSR to ELL format.
// rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); // A sparse matrix // 1 0 3 4 // 0 0 5 1 // 0 2 0 0 // 4 0 0 8 rocsparse_int hAptr[5] = {0, 3, 5, 6, 8}; rocsparse_int hAcol[8] = {0, 2, 3, 2, 3, 1, 0, 3}; double hAval[8] = {1.0, 3.0, 4.0, 5.0, 1.0, 2.0, 4.0, 8.0}; rocsparse_int m = 4; rocsparse_int n = 4; rocsparse_int nnz = 8; double halpha = 1.0; double hbeta = 0.0; double hx[4] = {1.0, 2.0, 3.0, 4.0}; double hy[4] = {4.0, 5.0, 6.0, 7.0}; // Matrix descriptors rocsparse_mat_descr descrA; rocsparse_create_mat_descr(&descrA); rocsparse_mat_descr descrB; rocsparse_create_mat_descr(&descrB); // Offload data to device rocsparse_int* dAptr = NULL; rocsparse_int* dAcol = NULL; double* dAval = NULL; double* dx = NULL; double* dy = NULL; hipMalloc((void**)&dAptr, sizeof(rocsparse_int) * (m + 1)); hipMalloc((void**)&dAcol, sizeof(rocsparse_int) * nnz); hipMalloc((void**)&dAval, sizeof(double) * nnz); hipMalloc((void**)&dx, sizeof(double) * n); hipMalloc((void**)&dy, sizeof(double) * m); hipMemcpy(dAptr, hAptr, sizeof(rocsparse_int) * (m + 1), hipMemcpyHostToDevice); hipMemcpy(dAcol, hAcol, sizeof(rocsparse_int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dAval, hAval, sizeof(double) * nnz, hipMemcpyHostToDevice); hipMemcpy(dx, hx, sizeof(double) * n, hipMemcpyHostToDevice); // Convert CSR matrix to ELL format rocsparse_int* dBcol = NULL; double* dBval = NULL; // Determine ELL width rocsparse_int ell_width; rocsparse_csr2ell_width(handle, m, descrA, dAptr, descrB, &ell_width); // Allocate memory for ELL storage format hipMalloc((void**)&dBcol, sizeof(rocsparse_int) * ell_width * m); hipMalloc((void**)&dBval, sizeof(double) * ell_width * m); // Convert matrix from CSR to ELL rocsparse_dcsr2ell(handle, m, descrA, dAval, dAptr, dAcol, descrB, ell_width, dBval, dBcol); // Clean up CSR structures hipFree(dAptr); hipFree(dAcol); hipFree(dAval); // Call rocsparse ellmv rocsparse_dellmv(handle, rocsparse_operation_none, m, n, &halpha, descrB, dBval, dBcol, ell_width, dx, &hbeta, dy); // Copy result back to host hipMemcpy(hy, dy, sizeof(double) * m, hipMemcpyDeviceToHost); // Clear up on device rocsparse_destroy_mat_descr(descrA); rocsparse_destroy_mat_descr(descrB); rocsparse_destroy_handle(handle); hipFree(dBcol); hipFree(dBval); hipFree(dx); hipFree(dy);
Note
This function does not produce deterministic results when A is transposed.
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
trans – [in] matrix operation type.
m – [in] number of rows of the sparse ELL matrix.
n – [in] number of columns of the sparse ELL matrix.
alpha – [in] scalar \(\alpha\).
descr – [in] descriptor of the sparse ELL matrix. Currently, only rocsparse_matrix_type_general is supported.
ell_val – [in] array that contains the elements of the sparse ELL matrix. Padded elements should be zero.
ell_col_ind – [in] array that contains the column indices of the sparse ELL matrix. Padded column indices should be -1.
ell_width – [in] number of non-zero elements per row of the sparse ELL matrix.
x – [in] array of
n
elements ( \(op(A) == A\)) orm
elements ( \(op(A) == A^T\) or \(op(A) == A^H\)).beta – [in] scalar \(\beta\).
y – [inout] array of
m
elements ( \(op(A) == A\)) orn
elements ( \(op(A) == A^T\) or \(op(A) == A^H\)).
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
m
,n
orell_width
is invalid.rocsparse_status_invalid_pointer –
descr
,alpha
,ell_val
,ell_col_ind
,x
,beta
ory
pointer is invalid.rocsparse_status_not_implemented – rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_hybmv()#
-
rocsparse_status rocsparse_shybmv(rocsparse_handle handle, rocsparse_operation trans, const float *alpha, const rocsparse_mat_descr descr, const rocsparse_hyb_mat hyb, const float *x, const float *beta, float *y)#
-
rocsparse_status rocsparse_dhybmv(rocsparse_handle handle, rocsparse_operation trans, const double *alpha, const rocsparse_mat_descr descr, const rocsparse_hyb_mat hyb, const double *x, const double *beta, double *y)#
-
rocsparse_status rocsparse_chybmv(rocsparse_handle handle, rocsparse_operation trans, const rocsparse_float_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_hyb_mat hyb, const rocsparse_float_complex *x, const rocsparse_float_complex *beta, rocsparse_float_complex *y)#
-
rocsparse_status rocsparse_zhybmv(rocsparse_handle handle, rocsparse_operation trans, const rocsparse_double_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_hyb_mat hyb, const rocsparse_double_complex *x, const rocsparse_double_complex *beta, rocsparse_double_complex *y)#
Sparse matrix vector multiplication using HYB storage format.
rocsparse_hybmv
multiplies the scalar \(\alpha\) with a sparse \(m \times n\) matrix, defined in HYB storage format, and the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == rocsparse_operation_none} \\ A^T, & \text{if trans == rocsparse_operation_transpose} \\ A^H, & \text{if trans == rocsparse_operation_conjugate_transpose} \end{array} \right. \end{split}\]- Example
This example performs a sparse matrix vector multiplication in HYB format. Also demonstrate conversion from CSR to HYB.
// rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); // A sparse matrix // 1 0 3 4 // 0 0 5 1 // 0 2 0 0 // 4 0 0 8 rocsparse_int hAptr[5] = {0, 3, 5, 6, 8}; rocsparse_int hAcol[8] = {0, 2, 3, 2, 3, 1, 0, 3}; double hAval[8] = {1.0, 3.0, 4.0, 5.0, 1.0, 2.0, 4.0, 8.0}; rocsparse_int m = 4; rocsparse_int n = 4; rocsparse_int nnz = 8; double halpha = 1.0; double hbeta = 0.0; double hx[4] = {1.0, 2.0, 3.0, 4.0}; double hy[4] = {4.0, 5.0, 6.0, 7.0}; // Matrix descriptor rocsparse_mat_descr descrA; rocsparse_create_mat_descr(&descrA); // Offload data to device rocsparse_int* dAptr = NULL; rocsparse_int* dAcol = NULL; double* dAval = NULL; double* dx = NULL; double* dy = NULL; hipMalloc((void**)&dAptr, sizeof(rocsparse_int) * (m + 1)); hipMalloc((void**)&dAcol, sizeof(rocsparse_int) * nnz); hipMalloc((void**)&dAval, sizeof(double) * nnz); hipMalloc((void**)&dx, sizeof(double) * n); hipMalloc((void**)&dy, sizeof(double) * m); hipMemcpy(dAptr, hAptr, sizeof(rocsparse_int) * (m + 1), hipMemcpyHostToDevice); hipMemcpy(dAcol, hAcol, sizeof(rocsparse_int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dAval, hAval, sizeof(double) * nnz, hipMemcpyHostToDevice); hipMemcpy(dx, hx, sizeof(double) * n, hipMemcpyHostToDevice); // Convert CSR matrix to HYB format rocsparse_hyb_mat hybA; rocsparse_create_hyb_mat(&hybA); rocsparse_dcsr2hyb(handle, m, n, descrA, dAval, dAptr, dAcol, hybA, 0, rocsparse_hyb_partition_auto); // Clean up CSR structures hipFree(dAptr); hipFree(dAcol); hipFree(dAval); // Call rocsparse hybmv rocsparse_dhybmv(handle, rocsparse_operation_none, &halpha, descrA, hybA, dx, &hbeta, dy); // Copy result back to host hipMemcpy(hy, dy, sizeof(double) * m, hipMemcpyDeviceToHost); // Clear up on device rocsparse_destroy_hyb_mat(hybA); rocsparse_destroy_mat_descr(descrA); rocsparse_destroy_handle(handle); hipFree(dx); hipFree(dy);
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
trans – [in] matrix operation type.
alpha – [in] scalar \(\alpha\).
descr – [in] descriptor of the sparse HYB matrix. Currently, only rocsparse_matrix_type_general is supported.
hyb – [in] matrix in HYB storage format.
x – [in] array of
n
elements ( \(op(A) == A\)) orm
elements ( \(op(A) == A^T\) or \(op(A) == A^H\)).beta – [in] scalar \(\beta\).
y – [inout] array of
m
elements ( \(op(A) == A\)) orn
elements ( \(op(A) == A^T\) or \(op(A) == A^H\)).
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
hyb
structure was not initialized with valid matrix sizes.rocsparse_status_invalid_pointer –
descr
,alpha
,hyb
,x
,beta
ory
pointer is invalid.rocsparse_status_invalid_value –
hyb
structure was not initialized with a valid partitioning type.rocsparse_status_arch_mismatch – the device is not supported.
rocsparse_status_memory_error – the buffer could not be allocated.
rocsparse_status_internal_error – an internal error occurred.
rocsparse_status_not_implemented –
trans
!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_gebsrmv()#
-
rocsparse_status rocsparse_sgebsrmv(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const float *alpha, const rocsparse_mat_descr descr, const float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int row_block_dim, rocsparse_int col_block_dim, const float *x, const float *beta, float *y)#
-
rocsparse_status rocsparse_dgebsrmv(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const double *alpha, const rocsparse_mat_descr descr, const double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int row_block_dim, rocsparse_int col_block_dim, const double *x, const double *beta, double *y)#
-
rocsparse_status rocsparse_cgebsrmv(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const rocsparse_float_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int row_block_dim, rocsparse_int col_block_dim, const rocsparse_float_complex *x, const rocsparse_float_complex *beta, rocsparse_float_complex *y)#
-
rocsparse_status rocsparse_zgebsrmv(rocsparse_handle handle, rocsparse_direction dir, rocsparse_operation trans, rocsparse_int mb, rocsparse_int nb, rocsparse_int nnzb, const rocsparse_double_complex *alpha, const rocsparse_mat_descr descr, const rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int row_block_dim, rocsparse_int col_block_dim, const rocsparse_double_complex *x, const rocsparse_double_complex *beta, rocsparse_double_complex *y)#
Sparse matrix vector multiplication using GEBSR storage format.
rocsparse_gebsrmv
multiplies the scalar \(\alpha\) with a sparse \((mb \cdot \text{row_block_dim}) \times (nb \cdot \text{col_block_dim})\) matrix, defined in GEBSR storage format, and the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == rocsparse_operation_none} \\ A^T, & \text{if trans == rocsparse_operation_transpose} \\ A^H, & \text{if trans == rocsparse_operation_conjugate_transpose} \end{array} \right. \end{split}\]- Example
This example performs a sparse matrix vector multiplication in GEBSR format.
// rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); // alpha * ( 1.0 0.0 2.0 ) * ( 1.0 ) + beta * ( 4.0 ) = ( 31.1 ) // ( 3.0 0.0 4.0 ) * ( 2.0 ) ( 5.0 ) = ( 62.0 ) // ( 5.0 6.0 0.0 ) * ( 3.0 ) ( 6.0 ) = ( 70.7 ) // ( 7.0 0.0 8.0 ) * ( 7.0 ) = ( 123.8 ) // GEBSR block dimensions rocsparse_int row_block_dim = 2; rocsparse_int col_block_dim = 3; // Number of block rows and columns rocsparse_int mb = 2; rocsparse_int nb = 1; // Number of non-zero blocks rocsparse_int nnzb = 2; // BSR row pointers rocsparse_int hbsr_row_ptr[3] = {0, 1, 2}; // BSR column indices rocsparse_int hbsr_col_ind[2] = {0, 0}; // BSR values double hbsr_val[16] = {1.0, 3.0, 0.0, 0.0, 2.0, 4.0, 5.0, 7.0, 6.0, 0.0, 0.0, 8.0}; // Block storage in column major rocsparse_direction dir = rocsparse_direction_column; // Transposition of the matrix rocsparse_operation trans = rocsparse_operation_none; // Scalar alpha and beta double alpha = 3.7; double beta = 1.3; // x and y double hx[4] = {1.0, 2.0, 3.0, 0.0}; double hy[4] = {4.0, 5.0, 6.0, 7.0}; // Matrix descriptor rocsparse_mat_descr descr; rocsparse_create_mat_descr(&descr); // Offload data to device rocsparse_int* dbsr_row_ptr; rocsparse_int* dbsr_col_ind; double* dbsr_val; double* dx; double* dy; hipMalloc((void**)&dbsr_row_ptr, sizeof(rocsparse_int) * (mb + 1)); hipMalloc((void**)&dbsr_col_ind, sizeof(rocsparse_int) * nnzb); hipMalloc((void**)&dbsr_val, sizeof(double) * nnzb * row_block_dim * col_block_dim); hipMalloc((void**)&dx, sizeof(double) * nb * col_block_dim); hipMalloc((void**)&dy, sizeof(double) * mb * row_block_dim); hipMemcpy(dbsr_row_ptr, hbsr_row_ptr, sizeof(rocsparse_int) * (mb + 1), hipMemcpyHostToDevice); hipMemcpy(dbsr_col_ind, hbsr_col_ind, sizeof(rocsparse_int) * nnzb, hipMemcpyHostToDevice); hipMemcpy(dbsr_val, hbsr_val, sizeof(double) * nnzb * row_block_dim * col_block_dim, hipMemcpyHostToDevice); hipMemcpy(dx, hx, sizeof(double) * nb * col_block_dim, hipMemcpyHostToDevice); hipMemcpy(dy, hy, sizeof(double) * mb * row_block_dim, hipMemcpyHostToDevice); // Call dbsrmv to perform y = alpha * A x + beta * y rocsparse_dgebsrmv(handle, dir, trans, mb, nb, nnzb, &alpha, descr, dbsr_val, dbsr_row_ptr, dbsr_col_ind, row_block_dim, col_block_dim, dx, &beta, dy); // Copy result back to host hipMemcpy(hy, dy, sizeof(double) * mb * row_block_dim, hipMemcpyDeviceToHost); // Clear rocSPARSE rocsparse_destroy_mat_descr(descr); rocsparse_destroy_handle(handle); // Clear device memory hipFree(dbsr_row_ptr); hipFree(dbsr_col_ind); hipFree(dbsr_val); hipFree(dx); hipFree(dy);
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
Currently, only
trans
== rocsparse_operation_none is supported.Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
dir – [in] matrix storage of GEBSR blocks.
trans – [in] matrix operation type.
mb – [in] number of block rows of the sparse GEBSR matrix.
nb – [in] number of block columns of the sparse GEBSR matrix.
nnzb – [in] number of non-zero blocks of the sparse GEBSR matrix.
alpha – [in] scalar \(\alpha\).
descr – [in] descriptor of the sparse GEBSR matrix. Currently, only rocsparse_matrix_type_general is supported.
bsr_val – [in] array of
nnzb
blocks of the sparse GEBSR matrix.bsr_row_ptr – [in] array of
mb+1
elements that point to the start of every block row of the sparse GEBSR matrix.bsr_col_ind – [in] array of
nnz
containing the block column indices of the sparse GEBSR matrix.row_block_dim – [in] row block dimension of the sparse GEBSR matrix.
col_block_dim – [in] column block dimension of the sparse GEBSR matrix.
x – [in] array of
nb*col_block_dim
elements ( \(op(A) = A\)) ormb*row_block_dim
elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).beta – [in] scalar \(\beta\).
y – [inout] array of
mb*row_block_dim
elements ( \(op(A) = A\)) ornb*col_block_dim
elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
mb
,nb
,nnzb
,row_block_dim
orcol_block_dim
is invalid.rocsparse_status_invalid_pointer –
descr
,alpha
,bsr_val
,bsr_row_ind
,bsr_col_ind
,x
,beta
ory
pointer is invalid.rocsparse_status_arch_mismatch – the device is not supported.
rocsparse_status_not_implemented –
trans
!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_gemvi_buffer_size()#
-
rocsparse_status rocsparse_sgemvi_buffer_size(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, size_t *buffer_size)#
-
rocsparse_status rocsparse_dgemvi_buffer_size(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, size_t *buffer_size)#
-
rocsparse_status rocsparse_cgemvi_buffer_size(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, size_t *buffer_size)#
-
rocsparse_status rocsparse_zgemvi_buffer_size(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, rocsparse_int nnz, size_t *buffer_size)#
Dense matrix sparse vector multiplication.
rocsparse_gemvi_buffer_size
returns the size of the temporary storage buffer required by rocsparse_sgemvi(), rocsparse_dgemvi(), rocsparse_cgemvi() or rocsparse_zgemvi(). The temporary storage buffer must be allocated by the user.Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
trans – [in] matrix operation type.
m – [in] number of rows of the dense matrix.
n – [in] number of columns of the dense matrix.
nnz – [in] number of non-zero entries in the sparse vector.
buffer_size – [out] temporary storage buffer size.
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
m
,n
, ornnz
is invalid.rocsparse_status_invalid_pointer –
buffer_size
pointer is invalid.rocsparse_status_not_implemented –
trans
!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
rocsparse_gemvi()#
-
rocsparse_status rocsparse_sgemvi(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, const float *alpha, const float *A, rocsparse_int lda, rocsparse_int nnz, const float *x_val, const rocsparse_int *x_ind, const float *beta, float *y, rocsparse_index_base idx_base, void *temp_buffer)#
-
rocsparse_status rocsparse_dgemvi(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, const double *alpha, const double *A, rocsparse_int lda, rocsparse_int nnz, const double *x_val, const rocsparse_int *x_ind, const double *beta, double *y, rocsparse_index_base idx_base, void *temp_buffer)#
-
rocsparse_status rocsparse_cgemvi(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, const rocsparse_float_complex *alpha, const rocsparse_float_complex *A, rocsparse_int lda, rocsparse_int nnz, const rocsparse_float_complex *x_val, const rocsparse_int *x_ind, const rocsparse_float_complex *beta, rocsparse_float_complex *y, rocsparse_index_base idx_base, void *temp_buffer)#
-
rocsparse_status rocsparse_zgemvi(rocsparse_handle handle, rocsparse_operation trans, rocsparse_int m, rocsparse_int n, const rocsparse_double_complex *alpha, const rocsparse_double_complex *A, rocsparse_int lda, rocsparse_int nnz, const rocsparse_double_complex *x_val, const rocsparse_int *x_ind, const rocsparse_double_complex *beta, rocsparse_double_complex *y, rocsparse_index_base idx_base, void *temp_buffer)#
Dense matrix sparse vector multiplication.
rocsparse_gemvi
multiplies the scalar \(\alpha\) with a dense \(m \times n\) matrix \(A\) and the sparse vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == rocsparse_operation_none} \\ A^T, & \text{if trans == rocsparse_operation_transpose} \\ A^H, & \text{if trans == rocsparse_operation_conjugate_transpose} \end{array} \right. \end{split}\]rocsparse_gemvi
requires a user allocated temporary buffer. Its size is returned by rocsparse_sgemvi_buffer_size(), rocsparse_dgemvi_buffer_size(), rocsparse_cgemvi_buffer_size() or rocsparse_zgemvi_buffer_size().- Example
// rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); rocsparse_mat_descr descr; rocsparse_create_mat_descr(&descr); // Number of rows and columns rocsparse_int m = 3; rocsparse_int n = 5; // Dense matrix A in column-major rocsparse_int lda = m; double hA[15] = {9, 14, 19, 10, 15, 20, 11, 16, 21, 12, 17, 22, 13, 18, 23}; // Number of non-zero entries rocsparse_int nnz = 3; // Vector x indices rocsparse_int hx_ind[3] = {0, 1, 3}; // Vector x values double hx_val[3] = {1, 2, 3}; // Vector y values double hy[3] = {4, 5, 6}; // Scalar alpha and beta double alpha = 3.7; double beta = 1.3; // Matrix operation rocsparse_operation trans = rocsparse_operation_none; // Index base rocsparse_index_base base = rocsparse_index_base_zero; // Offload data to device double* dA; rocsparse_int* dx_ind; double* dx_val; double* dy; hipMalloc((void**)&dA, sizeof(double) * m * n); hipMalloc((void**)&dx_ind, sizeof(rocsparse_int) * nnz); hipMalloc((void**)&dx_val, sizeof(double) * nnz); hipMalloc((void**)&dy, sizeof(double) * m); hipMemcpy(dA, hA, sizeof(double) * m * n, hipMemcpyHostToDevice); hipMemcpy(dx_ind, hx_ind, sizeof(rocsparse_int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dx_val, hx_val, sizeof(double) * nnz, hipMemcpyHostToDevice); hipMemcpy(dy, hy, sizeof(double) * m, hipMemcpyHostToDevice); // Obtain buffer size size_t buffer_size; rocsparse_dgemvi_buffer_size(handle, trans, m, n, nnz, &buffer_size); void* buffer; hipMalloc(&buffer, buffer_size); // Call dgemvi rocsparse_dgemvi(handle, trans, m, n, &alpha, dA, lda, nnz, dx_val, dx_ind, &beta, dy, base, buffer); // Copy result back to host hipMemcpy(hy, dy, sizeof(double) * m, hipMemcpyDeviceToHost); // Clear rocSPARSE rocsparse_destroy_handle(handle); // Clear device memory hipFree(dA); hipFree(dx_ind); hipFree(dx_val); hipFree(dy); hipFree(buffer);
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
Currently, only
trans
== rocsparse_operation_none is supported.Note
This routine supports execution in a hipGraph context.
- Parameters:
handle – [in] handle to the rocsparse library context queue.
trans – [in] matrix operation type.
m – [in] number of rows of the dense matrix.
n – [in] number of columns of the dense matrix.
alpha – [in] scalar \(\alpha\).
A – [in] pointer to the dense matrix.
lda – [in] leading dimension of the dense matrix
nnz – [in] number of non-zero entries in the sparse vector
x_val – [in] array of
nnz
elements containing the values of the sparse vectorx_ind – [in] array of
nnz
elements containing the indices of the sparse vectorbeta – [in] scalar \(\beta\).
y – [inout] array of
m
elements ( \(op(A) == A\)) orn
elements ( \(op(A) == A^T\) or \(op(A) == A^H\)).idx_base – [in] rocsparse_index_base_zero or rocsparse_index_base_one.
temp_buffer – [in] temporary storage buffer
- Return values:
rocsparse_status_success – the operation completed successfully.
rocsparse_status_invalid_handle – the library context was not initialized.
rocsparse_status_invalid_size –
m
,n
,lda
ornnz
is invalid.rocsparse_status_invalid_pointer –
alpha
,A
,x_val
,x_ind
,beta
,y
ortemp_buffer
pointer is invalid.rocsparse_status_not_implemented –
trans
!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.