Sparse level 2 functions

Contents

Sparse level 2 functions#

This module contains all sparse level 2 routines.

The sparse level 2 routines describe operations between a matrix in sparse format and a vector in dense format.

hipsparseXcsrmv()#

hipsparseStatus_t hipsparseScsrmv(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, const float *alpha, const hipsparseMatDescr_t descrA, const float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const float *x, const float *beta, float *y)#
hipsparseStatus_t hipsparseDcsrmv(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, const double *alpha, const hipsparseMatDescr_t descrA, const double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const double *x, const double *beta, double *y)#
hipsparseStatus_t hipsparseCcsrmv(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, const hipComplex *alpha, const hipsparseMatDescr_t descrA, const hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const hipComplex *x, const hipComplex *beta, hipComplex *y)#
hipsparseStatus_t hipsparseZcsrmv(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, const hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const hipDoubleComplex *x, const hipDoubleComplex *beta, hipDoubleComplex *y)#

Sparse matrix vector multiplication using CSR storage format.

hipsparseXcsrmv 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 transA == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ A^T, & \text{if transA == HIPSPARSE_OPERATION_TRANSPOSE} \\ A^H, & \text{if transA == HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE} \end{array} \right. \end{split}\]

for(i = 0; i < m; ++i)
{
    y[i] = beta * y[i];

    for(j = csrRowPtr[i]; j < csrRowPtr[i + 1]; ++j)
    {
        y[i] = y[i] + alpha * csrVal[j] * x[csrColInd[j]];
    }
}

Deprecated:

This function is deprecated when using the CUDA backend (CUDA 10.0+) and will be removed in CUDA 11.0. This deprecation does not apply to the ROCm backend.

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 transA == HIPSPARSE_OPERATION_NON_TRANSPOSE is supported.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • transA[in] matrix operation type.

  • m[in] number of rows of the sparse CSR matrix. Must be non-negative.

  • n[in] number of columns of the sparse CSR matrix. Must be non-negative.

  • nnz[in] number of non-zero entries of the sparse CSR matrix. Must be non-negative.

  • alpha[in] scalar \(\alpha\).

  • descrA[in] descriptor of the sparse CSR matrix. Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • csrSortedValA[in] array of nnz elements of the sparse CSR matrix.

  • csrSortedRowPtrA[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix.

  • csrSortedColIndA[in] array of nnz elements containing the column indices of the sparse CSR matrix.

  • x[in] array of n elements ( \(op(A) == A\)) or m elements ( \(op(A) == A^T\) or \(op(A) == A^H\)).

  • beta[in] scalar \(\beta\).

  • y[inout] array of m elements ( \(op(A) == A\)) or n elements ( \(op(A) == A^T\) or \(op(A) == A^H\)).

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, descrA, alpha or beta is nullptr, m, n or nnz is negative, or csrSortedValA, csrSortedRowPtrA, csrSortedColIndA, x or y is nullptr.

  • HIPSPARSE_STATUS_ARCH_MISMATCH – the device is not supported.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t is not HIPSPARSE_MATRIX_TYPE_GENERAL.

 1int main(int argc, char* argv[])
 2{
 3    // hipSPARSE handle
 4    hipsparseHandle_t handle;
 5    HIPSPARSE_CHECK(hipsparseCreate(&handle));
 6
 7    // alpha * ( 1.0  0.0  2.0 ) * ( 1.0 ) + beta * ( 4.0 ) = (  31.1 )
 8    //         ( 3.0  0.0  4.0 ) * ( 2.0 )          ( 5.0 ) = (  62.0 )
 9    //         ( 5.0  6.0  0.0 ) * ( 3.0 )          ( 6.0 ) = (  70.7 )
10    //         ( 7.0  0.0  8.0 ) *                  ( 7.0 ) = ( 123.8 )
11
12    const int m   = 4;
13    const int n   = 3;
14    const int nnz = 8;
15
16    // CSR row pointers
17    std::vector<int> hcsrRowPtr = {0, 2, 4, 6, 8};
18
19    // CSR column indices
20    std::vector<int> hcsrColInd = {0, 2, 0, 2, 0, 1, 0, 2};
21
22    // CSR values
23    std::vector<double> hcsrVal = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0};
24
25    // Transposition of the matrix
26    hipsparseOperation_t trans = HIPSPARSE_OPERATION_NON_TRANSPOSE;
27
28    // Scalar alpha and beta
29    double alpha = 3.7;
30    double beta  = 1.3;
31
32    // x and y
33    std::vector<double> hx = {1.0, 2.0, 3.0};
34    std::vector<double> hy = {4.0, 5.0, 6.0, 7.0};
35
36    // Matrix descriptor
37    hipsparseMatDescr_t descr;
38    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descr));
39
40    // Offload data to device
41    int*    dcsrRowPtr;
42    int*    dcsrColInd;
43    double* dcsrVal;
44    double* dx;
45    double* dy;
46
47    HIP_CHECK(hipMalloc((void**)&dcsrRowPtr, sizeof(int) * (m + 1)));
48    HIP_CHECK(hipMalloc((void**)&dcsrColInd, sizeof(int) * nnz));
49    HIP_CHECK(hipMalloc((void**)&dcsrVal, sizeof(double) * nnz));
50    HIP_CHECK(hipMalloc((void**)&dx, sizeof(double) * n));
51    HIP_CHECK(hipMalloc((void**)&dy, sizeof(double) * m));
52
53    HIP_CHECK(
54        hipMemcpy(dcsrRowPtr, hcsrRowPtr.data(), sizeof(int) * (m + 1), hipMemcpyHostToDevice));
55    HIP_CHECK(hipMemcpy(dcsrColInd, hcsrColInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
56    HIP_CHECK(hipMemcpy(dcsrVal, hcsrVal.data(), sizeof(double) * nnz, hipMemcpyHostToDevice));
57    HIP_CHECK(hipMemcpy(dx, hx.data(), sizeof(double) * n, hipMemcpyHostToDevice));
58    HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(double) * m, hipMemcpyHostToDevice));
59
60    // Call dcsrmv to perform y = alpha * A x + beta * y
61    HIPSPARSE_CHECK(hipsparseDcsrmv(
62        handle, trans, m, n, nnz, &alpha, descr, dcsrVal, dcsrRowPtr, dcsrColInd, dx, &beta, dy));
63
64    // Copy result back to host
65    HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(double) * m, hipMemcpyDeviceToHost));
66
67    std::cout << "hy" << std::endl;
68    for(int i = 0; i < m; i++)
69    {
70        std::cout << hy[i] << " ";
71    }
72    std::cout << std::endl;
73
74    // Clear hipSPARSE
75    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descr));
76    HIPSPARSE_CHECK(hipsparseDestroy(handle));
77
78    // Clear device memory
79    HIP_CHECK(hipFree(dcsrRowPtr));
80    HIP_CHECK(hipFree(dcsrColInd));
81    HIP_CHECK(hipFree(dcsrVal));
82    HIP_CHECK(hipFree(dx));
83    HIP_CHECK(hipFree(dy));
84
85    return 0;
86}

hipsparseXcsrsv2_zeroPivot()#

hipsparseStatus_t hipsparseXcsrsv2_zeroPivot(hipsparseHandle_t handle, csrsv2Info_t info, int *position)#

hipsparseXcsrsv2_zeroPivot returns HIPSPARSE_STATUS_ZERO_PIVOT, if either a structural or numerical zero has been found during hipsparseScsrsv2_solve(), hipsparseDcsrsv2_solve(), hipsparseCcsrsv2_solve() or hipsparseZcsrsv2_solve() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored in position, 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 HIPSPARSE_STATUS_SUCCESS is returned instead.

Deprecated:

This function is deprecated when using the CUDA backend (CUDA 11.0+) and will be removed in CUDA 12.0. This deprecation does not apply to the ROCm backend.

Note

hipsparseXcsrsv2_zeroPivot is a blocking function. It might influence performance negatively.

Parameters:
  • handle[in] handle to the hipsparse 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:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, info or position is nullptr.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_ZERO_PIVOT – zero pivot has been found.

hipsparseXcsrsv2_bufferSize()#

hipsparseStatus_t hipsparseScsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDcsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCcsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZcsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, int *pBufferSizeInBytes)#

hipsparseXcsrsv2_bufferSize returns the size of the temporary storage buffer in bytes that is required by hipsparseScsrsv2_analysis() and hipsparseXcsrsv2_solve(). The temporary storage buffer must be allocated by the user.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • transA[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.

  • descrA[in] descriptor of the sparse CSR matrix.

  • csrSortedValA[in] array of nnz elements of the sparse CSR matrix.

  • csrSortedRowPtrA[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix.

  • csrSortedColIndA[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.

  • pBufferSizeInBytes[out] number of bytes of the temporary storage buffer required by hipsparseXcsrsv2_analysis() and hipsparseXcsrsv2_solve().

Return values:

hipsparseXcsrsv2_bufferSizeExt()#

hipsparseStatus_t hipsparseScsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDcsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCcsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZcsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, size_t *pBufferSizeInBytes)#

hipsparseXcsrsv2_bufferSizeExt returns the size of the temporary storage buffer in bytes that is required by hipsparseXcsrsv2_analysis() and hipsparseScsrsv2_solve(). The temporary storage buffer must be allocated by the user.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • transA[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.

  • descrA[in] descriptor of the sparse CSR matrix.

  • csrSortedValA[in] array of nnz elements of the sparse CSR matrix.

  • csrSortedRowPtrA[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix.

  • csrSortedColIndA[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.

  • pBufferSizeInBytes[out] number of bytes of the temporary storage buffer required by hipsparseXcsrsv2_analysis() and hipsparseXcsrsv2_solve().

Return values:

hipsparseXcsrsv2_analysis()#

hipsparseStatus_t hipsparseScsrsv2_analysis(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, const float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseDcsrsv2_analysis(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, const double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseCcsrsv2_analysis(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, const hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseZcsrsv2_analysis(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, const hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#

hipsparseXcsrsv2_analysis performs the analysis step for hipsparseXcsrsv2_solve(). It is expected that this function will be executed only once for a given matrix and particular operation type.

Note

This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • transA[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.

  • descrA[in] descriptor of the sparse CSR matrix.

  • csrSortedValA[in] array of nnz elements of the sparse CSR matrix.

  • csrSortedRowPtrA[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix.

  • csrSortedColIndA[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.

  • policy[in] HIPSPARSE_SOLVE_POLICY_NO_LEVEL or HIPSPARSE_SOLVE_POLICY_USE_LEVEL.

  • pBuffer[in] temporary storage buffer allocated by the user.

Return values:

hipsparseXcsrsv2_solve()#

hipsparseStatus_t hipsparseScsrsv2_solve(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const float *alpha, const hipsparseMatDescr_t descrA, const float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, const float *f, float *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseDcsrsv2_solve(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const double *alpha, const hipsparseMatDescr_t descrA, const double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, const double *f, double *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseCcsrsv2_solve(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipComplex *alpha, const hipsparseMatDescr_t descrA, const hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, const hipComplex *f, hipComplex *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseZcsrsv2_solve(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, const hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, const hipDoubleComplex *f, hipDoubleComplex *x, hipsparseSolvePolicy_t policy, void *pBuffer)#

Sparse triangular solve using CSR storage format.

hipsparseXcsrsv2_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 == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ A^T, & \text{if trans == HIPSPARSE_OPERATION_TRANSPOSE} \\ A^H, & \text{if trans == HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE} \end{array} \right. \end{split}\]

Performing the above operation requires three steps. First, the user calls hipsparseXcsrsv2_bufferSize() (or hipsparseXcsrsv2_bufferSizeExt()) which will determine the size of the required temporary storage buffer. The user then allocates this buffer and calls hipsparseXcsrsv2_analysis() which will perform analysis on the sparse matrix \(op(A)\). Finally, the user completes the computation by calling hipsparseXcsrsv2_solve. The buffer size, buffer allocation, and analysis only need to be called once for a given sparse matrix \(op(A)\) while the computation stage can be repeatedly used with different \(x\) and \(y\) vectors. Once all calls to hipsparseXcsrsv2_solve are complete, the temporary buffer can be deallocated.

Solving a triangular system involves division by the diagonal elements. This means that if the sparse matrix is missing the diagonal entry (referred to as a structural zero) or the diagonal entry is zero (referred to as a numerical zero) then a division by zero would occur. hipsparseXcsrsv2_solve tracks the location of the first zero pivot (either numerical or structural zero). The zero pivot status can be checked calling hipsparseXcsrsv2_zeroPivot(). If hipsparseXcsrsv2_zeroPivot() returns HIPSPARSE_STATUS_SUCCESS, then no zero pivot was found and therefore the matrix does not have a structural or numerical zero.

The user can specify that the sparse matrix should be interpreted as having ones on the diagonal by setting the diagonal type on the descriptor descrA to HIPSPARSE_DIAG_TYPE_UNIT using hipsparseSetMatDiagType. If hipsparseDiagType_t == HIPSPARSE_DIAG_TYPE_UNIT, no zero pivot will be reported, even if \(A_{j,j} = 0\) for some \(j\).

The sparse CSR matrix passed to hipsparseXcsrsv2_solve does not actually have to be a triangular matrix. Instead the triangular upper or lower part of the sparse matrix is solved based on hipsparseFillMode_t set on the descriptor descrA. If the fill mode is set to HIPSPARSE_FILL_MODE_LOWER, then the lower triangular matrix is solved. If the fill mode is set to HIPSPARSE_FILL_MODE_UPPER then the upper triangular matrix is solved.

Note

The sparse CSR matrix has to be sorted. This can be achieved by calling hipsparseXcsrsort().

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 transA == HIPSPARSE_OPERATION_NON_TRANSPOSE and transA == HIPSPARSE_OPERATION_TRANSPOSE is supported.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • transA[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\).

  • descrA[in] descriptor of the sparse CSR matrix.

  • csrSortedValA[in] array of nnz elements of the sparse CSR matrix.

  • csrSortedRowPtrA[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix.

  • csrSortedColIndA[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.

  • f[in] array of m elements, holding the right-hand side.

  • x[out] array of m elements, holding the solution.

  • policy[in] HIPSPARSE_SOLVE_POLICY_NO_LEVEL or HIPSPARSE_SOLVE_POLICY_USE_LEVEL.

  • pBuffer[in] temporary storage buffer allocated by the user.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, nnz, descrA, alpha, csrSortedValA, csrSortedRowPtrA, csrSortedColIndA, f or x is invalid.

  • HIPSPARSE_STATUS_ARCH_MISMATCH – the device is not supported.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDtransA == HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE or hipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

  1int main(int argc, char* argv[])
  2{
  3    // hipSPARSE handle
  4    hipsparseHandle_t handle;
  5    HIPSPARSE_CHECK(hipsparseCreate(&handle));
  6
  7    // alpha * ( 1.0  0.0  2.0  0.0 ) * ( x_0 ) = ( 32.0 )
  8    //         ( 3.0  2.0  4.0  1.0 ) * ( x_1 ) = ( 14.7 )
  9    //         ( 5.0  6.0  1.0  3.0 ) * ( x_2 ) = ( 33.6 )
 10    //         ( 7.0  0.0  8.0  0.6 ) * ( x_3 ) = ( 10.0 )
 11
 12    const int m   = 4;
 13    const int nnz = 13;
 14
 15    // CSR row pointers
 16    std::vector<int> hcsrRowPtr = {0, 2, 6, 10, 13};
 17
 18    // CSR column indices
 19    std::vector<int> hcsrColInd = {0, 2, 0, 1, 2, 3, 0, 1, 2, 3, 0, 2, 3};
 20
 21    // CSR values
 22    std::vector<double> hcsrVal = {1.0, 2.0, 3.0, 2.0, 4.0, 1.0, 5.0, 6.0, 1.0, 3.0, 7.0, 8.0, 0.6};
 23
 24    // Transposition of the matrix
 25    hipsparseOperation_t   trans  = HIPSPARSE_OPERATION_NON_TRANSPOSE;
 26    hipsparseSolvePolicy_t policy = HIPSPARSE_SOLVE_POLICY_USE_LEVEL;
 27
 28    // Scalar alpha
 29    double alpha = 1.0;
 30
 31    // f and x
 32    std::vector<double> hf = {32.0, 14.7, 33.6, 10.0};
 33    std::vector<double> hx(m);
 34
 35    // Matrix descriptor
 36    hipsparseMatDescr_t descr;
 37    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descr));
 38
 39    // Set index base on descriptor
 40    HIPSPARSE_CHECK(hipsparseSetMatIndexBase(descr, HIPSPARSE_INDEX_BASE_ZERO));
 41
 42    // Set fill mode on descriptor
 43    HIPSPARSE_CHECK(hipsparseSetMatFillMode(descr, HIPSPARSE_FILL_MODE_LOWER));
 44
 45    // Set diag type on descriptor
 46    HIPSPARSE_CHECK(hipsparseSetMatDiagType(descr, HIPSPARSE_DIAG_TYPE_UNIT));
 47
 48    // Csrsv info
 49    csrsv2Info_t info;
 50    HIPSPARSE_CHECK(hipsparseCreateCsrsv2Info(&info));
 51
 52    // Offload data to device
 53    int*    dcsrRowPtr;
 54    int*    dcsrColInd;
 55    double* dcsrVal;
 56    double* df;
 57    double* dx;
 58
 59    HIP_CHECK(hipMalloc((void**)&dcsrRowPtr, sizeof(int) * (m + 1)));
 60    HIP_CHECK(hipMalloc((void**)&dcsrColInd, sizeof(int) * nnz));
 61    HIP_CHECK(hipMalloc((void**)&dcsrVal, sizeof(double) * nnz));
 62    HIP_CHECK(hipMalloc((void**)&df, sizeof(double) * m));
 63    HIP_CHECK(hipMalloc((void**)&dx, sizeof(double) * m));
 64
 65    HIP_CHECK(
 66        hipMemcpy(dcsrRowPtr, hcsrRowPtr.data(), sizeof(int) * (m + 1), hipMemcpyHostToDevice));
 67    HIP_CHECK(hipMemcpy(dcsrColInd, hcsrColInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
 68    HIP_CHECK(hipMemcpy(dcsrVal, hcsrVal.data(), sizeof(double) * nnz, hipMemcpyHostToDevice));
 69    HIP_CHECK(hipMemcpy(df, hf.data(), sizeof(double) * m, hipMemcpyHostToDevice));
 70
 71    int bufferSize = 0;
 72    HIPSPARSE_CHECK(hipsparseDcsrsv2_bufferSize(
 73        handle, trans, m, nnz, descr, dcsrVal, dcsrRowPtr, dcsrColInd, info, &bufferSize));
 74
 75    void* dbuffer = nullptr;
 76    HIP_CHECK(hipMalloc((void**)&dbuffer, bufferSize));
 77
 78    HIPSPARSE_CHECK(hipsparseDcsrsv2_analysis(
 79        handle, trans, m, nnz, descr, dcsrVal, dcsrRowPtr, dcsrColInd, info, policy, dbuffer));
 80
 81    // Call dcsrsv to perform alpha * A * x = f
 82    HIPSPARSE_CHECK(hipsparseDcsrsv2_solve(handle,
 83                                           trans,
 84                                           m,
 85                                           nnz,
 86                                           &alpha,
 87                                           descr,
 88                                           dcsrVal,
 89                                           dcsrRowPtr,
 90                                           dcsrColInd,
 91                                           info,
 92                                           df,
 93                                           dx,
 94                                           policy,
 95                                           dbuffer));
 96
 97    // Copy result back to host
 98    HIP_CHECK(hipMemcpy(hx.data(), dx, sizeof(double) * m, hipMemcpyDeviceToHost));
 99
100    std::cout << "hx" << std::endl;
101    for(int i = 0; i < m; i++)
102    {
103        std::cout << hx[i] << " ";
104    }
105    std::cout << std::endl;
106
107    // Clear hipSPARSE
108    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descr));
109    HIPSPARSE_CHECK(hipsparseDestroyCsrsv2Info(info));
110    HIPSPARSE_CHECK(hipsparseDestroy(handle));
111
112    // Clear device memory
113    HIP_CHECK(hipFree(dcsrRowPtr));
114    HIP_CHECK(hipFree(dcsrColInd));
115    HIP_CHECK(hipFree(dcsrVal));
116    HIP_CHECK(hipFree(df));
117    HIP_CHECK(hipFree(dx));
118    HIP_CHECK(hipFree(dbuffer));
119
120    return 0;
121}

hipsparseXhybmv()#

hipsparseStatus_t hipsparseShybmv(hipsparseHandle_t handle, hipsparseOperation_t transA, const float *alpha, const hipsparseMatDescr_t descrA, const hipsparseHybMat_t hybA, const float *x, const float *beta, float *y)#
hipsparseStatus_t hipsparseDhybmv(hipsparseHandle_t handle, hipsparseOperation_t transA, const double *alpha, const hipsparseMatDescr_t descrA, const hipsparseHybMat_t hybA, const double *x, const double *beta, double *y)#
hipsparseStatus_t hipsparseChybmv(hipsparseHandle_t handle, hipsparseOperation_t transA, const hipComplex *alpha, const hipsparseMatDescr_t descrA, const hipsparseHybMat_t hybA, const hipComplex *x, const hipComplex *beta, hipComplex *y)#
hipsparseStatus_t hipsparseZhybmv(hipsparseHandle_t handle, hipsparseOperation_t transA, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, const hipsparseHybMat_t hybA, const hipDoubleComplex *x, const hipDoubleComplex *beta, hipDoubleComplex *y)#

Sparse matrix vector multiplication using HYB storage format.

hipsparseXhybmv 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
\[ op(A) = \left\{ \begin{array}{ll} A, & \text{if transA == HIPSPARSE_OPERATION_NON_TRANSPOSE} \end{array} \right. \]

Deprecated:

This function is deprecated when using the CUDA backend (CUDA 10.0+) and will be removed in CUDA 11.0. This deprecation does not apply to the ROCm backend.

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 transA == HIPSPARSE_OPERATION_NON_TRANSPOSE is supported.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • transA[in] matrix operation type.

  • alpha[in] scalar \(\alpha\).

  • descrA[in] descriptor of the sparse HYB matrix. Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • hybA[in] matrix in HYB storage format.

  • x[in] array of n elements ( \(op(A) == A\)) or m elements ( \(op(A) == A^T\) or \(op(A) == A^H\)).

  • beta[in] scalar \(\beta\).

  • y[inout] array of m elements ( \(op(A) == A\)) or n elements ( \(op(A) == A^T\) or \(op(A) == A^H\)).

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, descrA, alpha, beta or hybA is nullptr, or x or y is nullptr.

  • HIPSPARSE_STATUS_ARCH_MISMATCH – the device is not supported.

  • HIPSPARSE_STATUS_ALLOC_FAILED – the buffer could not be allocated.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDtransA is not HIPSPARSE_OPERATION_NON_TRANSPOSE, or hipsparseMatrixType_t is not HIPSPARSE_MATRIX_TYPE_GENERAL.

 1int main(int argc, char* argv[])
 2{
 3    // hipSPARSE handle
 4    hipsparseHandle_t handle;
 5    HIPSPARSE_CHECK(hipsparseCreate(&handle));
 6
 7    // A sparse matrix
 8    // 1 0 3 4
 9    // 0 0 5 1
10    // 0 2 0 0
11    // 4 0 0 8
12    std::vector<int>    hAptr = {0, 3, 5, 6, 8};
13    std::vector<int>    hAcol = {0, 2, 3, 2, 3, 1, 0, 3};
14    std::vector<double> hAval = {1.0, 3.0, 4.0, 5.0, 1.0, 2.0, 4.0, 8.0};
15
16    int m   = 4;
17    int n   = 4;
18    int nnz = 8;
19
20    double halpha = 1.0;
21    double hbeta  = 0.0;
22
23    std::vector<double> hx = {1.0, 2.0, 3.0, 4.0};
24    std::vector<double> hy = {4.0, 5.0, 6.0, 7.0};
25
26    // Matrix descriptor
27    hipsparseMatDescr_t descrA;
28    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrA));
29
30    // Offload data to device
31    int*    dAptr = NULL;
32    int*    dAcol = NULL;
33    double* dAval = NULL;
34    double* dx    = NULL;
35    double* dy    = NULL;
36
37    HIP_CHECK(hipMalloc((void**)&dAptr, sizeof(int) * (m + 1)));
38    HIP_CHECK(hipMalloc((void**)&dAcol, sizeof(int) * nnz));
39    HIP_CHECK(hipMalloc((void**)&dAval, sizeof(double) * nnz));
40    HIP_CHECK(hipMalloc((void**)&dx, sizeof(double) * n));
41    HIP_CHECK(hipMalloc((void**)&dy, sizeof(double) * m));
42
43    HIP_CHECK(hipMemcpy(dAptr, hAptr.data(), sizeof(int) * (m + 1), hipMemcpyHostToDevice));
44    HIP_CHECK(hipMemcpy(dAcol, hAcol.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
45    HIP_CHECK(hipMemcpy(dAval, hAval.data(), sizeof(double) * nnz, hipMemcpyHostToDevice));
46    HIP_CHECK(hipMemcpy(dx, hx.data(), sizeof(double) * n, hipMemcpyHostToDevice));
47
48    // Convert CSR matrix to HYB format
49    hipsparseHybMat_t hybA;
50    HIPSPARSE_CHECK(hipsparseCreateHybMat(&hybA));
51
52    HIPSPARSE_CHECK(hipsparseDcsr2hyb(
53        handle, m, n, descrA, dAval, dAptr, dAcol, hybA, 0, HIPSPARSE_HYB_PARTITION_AUTO));
54
55    // Clean up CSR structures
56    HIP_CHECK(hipFree(dAptr));
57    HIP_CHECK(hipFree(dAcol));
58    HIP_CHECK(hipFree(dAval));
59
60    // Call hipsparse hybmv
61    HIPSPARSE_CHECK(hipsparseDhybmv(
62        handle, HIPSPARSE_OPERATION_NON_TRANSPOSE, &halpha, descrA, hybA, dx, &hbeta, dy));
63
64    // Copy result back to host
65    HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(double) * m, hipMemcpyDeviceToHost));
66
67    std::cout << "hy" << std::endl;
68    for(int i = 0; i < m; i++)
69    {
70        std::cout << hy[i] << " ";
71    }
72    std::cout << std::endl;
73
74    // Clear up on device
75    HIPSPARSE_CHECK(hipsparseDestroyHybMat(hybA));
76    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrA));
77    HIPSPARSE_CHECK(hipsparseDestroy(handle));
78
79    HIP_CHECK(hipFree(dx));
80    HIP_CHECK(hipFree(dy));
81
82    return 0;
83}

hipsparseXbsrmv()#

hipsparseStatus_t hipsparseSbsrmv(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nb, int nnzb, const float *alpha, const hipsparseMatDescr_t descrA, const float *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, const float *x, const float *beta, float *y)#
hipsparseStatus_t hipsparseDbsrmv(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nb, int nnzb, const double *alpha, const hipsparseMatDescr_t descrA, const double *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, const double *x, const double *beta, double *y)#
hipsparseStatus_t hipsparseCbsrmv(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nb, int nnzb, const hipComplex *alpha, const hipsparseMatDescr_t descrA, const hipComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, const hipComplex *x, const hipComplex *beta, hipComplex *y)#
hipsparseStatus_t hipsparseZbsrmv(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nb, int nnzb, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, const hipDoubleComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, const hipDoubleComplex *x, const hipDoubleComplex *beta, hipDoubleComplex *y)#

Sparse matrix vector multiplication using BSR storage format.

hipsparseXbsrmv multiplies the scalar \(\alpha\) with a sparse \(m \times n\) 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
\[ op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == HIPSPARSE_OPERATION_NON_TRANSPOSE} \end{array} \right. \]
and where \(m = mb \times blockDim\) and \(n= nb \times blockDim\).

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 transA == HIPSPARSE_OPERATION_NON_TRANSPOSE is supported.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • dirA[in] matrix storage of BSR blocks.

  • transA[in] matrix operation type.

  • mb[in] number of block rows of the sparse BSR matrix. Must be non-negative.

  • nb[in] number of block columns of the sparse BSR matrix. Must be non-negative.

  • nnzb[in] number of non-zero blocks of the sparse BSR matrix. Must be non-negative.

  • alpha[in] scalar \(\alpha\).

  • descrA[in] descriptor of the sparse BSR matrix. Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • bsrSortedValA[in] array of nnzb blocks of the sparse BSR matrix.

  • bsrSortedRowPtrA[in] array of mb+1 elements that point to the start of every block row of the sparse BSR matrix.

  • bsrSortedColIndA[in] array of nnzb elements containing the block column indices of the sparse BSR matrix.

  • blockDim[in] block dimension of the sparse BSR matrix. Must be positive.

  • x[in] array of nb*blockDim elements ( \(op(A) = A\)) or mb*blockDim elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).

  • beta[in] scalar \(\beta\).

  • y[inout] array of mb*blockDim elements ( \(op(A) = A\)) or nb*blockDim elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, descrA, alpha or beta is nullptr, mb, nb or nnzb is negative, blockDim is less than or equal to zero, or bsrSortedValA, bsrSortedRowPtrA, bsrSortedColIndA, x or y is nullptr.

  • HIPSPARSE_STATUS_ARCH_MISMATCH – the device is not supported.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDtransA is not HIPSPARSE_OPERATION_NON_TRANSPOSE, or hipsparseMatrixType_t is not HIPSPARSE_MATRIX_TYPE_GENERAL.

  1int main(int argc, char* argv[])
  2{
  3    // hipSPARSE handle
  4    hipsparseHandle_t handle;
  5    HIPSPARSE_CHECK(hipsparseCreate(&handle));
  6
  7    // alpha * ( 1.0  0.0  2.0 ) * ( 1.0 ) + beta * ( 4.0 ) = (  31.1 )
  8    //         ( 3.0  0.0  4.0 ) * ( 2.0 )          ( 5.0 ) = (  62.0 )
  9    //         ( 5.0  6.0  0.0 ) * ( 3.0 )          ( 6.0 ) = (  70.7 )
 10    //         ( 7.0  0.0  8.0 ) *                  ( 7.0 ) = ( 123.8 )
 11
 12    // BSR block dimension
 13    const int bsr_dim = 2;
 14
 15    // Number of block rows and columns
 16    const int mb = 2;
 17    const int nb = 2;
 18
 19    // Number of non-zero blocks
 20    const int nnzb = 4;
 21
 22    // Number of rows and columns
 23    const int m = mb * bsr_dim;
 24    const int n = nb * bsr_dim;
 25
 26    // BSR row pointers
 27    std::vector<int> hbsrRowPtr = {0, 2, 4};
 28
 29    // BSR column indices
 30    std::vector<int> hbsrColInd = {0, 1, 0, 1};
 31
 32    // BSR values
 33    std::vector<double> hbsrVal
 34        = {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};
 35
 36    // Block storage in column major
 37    hipsparseDirection_t dir = HIPSPARSE_DIRECTION_COLUMN;
 38
 39    // Transposition of the matrix
 40    hipsparseOperation_t trans = HIPSPARSE_OPERATION_NON_TRANSPOSE;
 41
 42    // Scalar alpha and beta
 43    double alpha = 3.7;
 44    double beta  = 1.3;
 45
 46    // x and y
 47    std::vector<double> hx = {1.0, 2.0, 3.0, 0.0};
 48    std::vector<double> hy = {4.0, 5.0, 6.0, 7.0};
 49
 50    // Matrix descriptor
 51    hipsparseMatDescr_t descr;
 52    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descr));
 53
 54    // Offload data to device
 55    int*    dbsrRowPtr;
 56    int*    dbsrColInd;
 57    double* dbsrVal;
 58    double* dx;
 59    double* dy;
 60
 61    HIP_CHECK(hipMalloc((void**)&dbsrRowPtr, sizeof(int) * (mb + 1)));
 62    HIP_CHECK(hipMalloc((void**)&dbsrColInd, sizeof(int) * nnzb));
 63    HIP_CHECK(hipMalloc((void**)&dbsrVal, sizeof(double) * nnzb * bsr_dim * bsr_dim));
 64    HIP_CHECK(hipMalloc((void**)&dx, sizeof(double) * n));
 65    HIP_CHECK(hipMalloc((void**)&dy, sizeof(double) * m));
 66
 67    HIP_CHECK(
 68        hipMemcpy(dbsrRowPtr, hbsrRowPtr.data(), sizeof(int) * (mb + 1), hipMemcpyHostToDevice));
 69    HIP_CHECK(hipMemcpy(dbsrColInd, hbsrColInd.data(), sizeof(int) * nnzb, hipMemcpyHostToDevice));
 70    HIP_CHECK(hipMemcpy(
 71        dbsrVal, hbsrVal.data(), sizeof(double) * nnzb * bsr_dim * bsr_dim, hipMemcpyHostToDevice));
 72    HIP_CHECK(hipMemcpy(dx, hx.data(), sizeof(double) * n, hipMemcpyHostToDevice));
 73    HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(double) * m, hipMemcpyHostToDevice));
 74
 75    // Call dbsrmv to perform y = alpha * A x + beta * y
 76    HIPSPARSE_CHECK(hipsparseDbsrmv(handle,
 77                                    dir,
 78                                    trans,
 79                                    mb,
 80                                    nb,
 81                                    nnzb,
 82                                    &alpha,
 83                                    descr,
 84                                    dbsrVal,
 85                                    dbsrRowPtr,
 86                                    dbsrColInd,
 87                                    bsr_dim,
 88                                    dx,
 89                                    &beta,
 90                                    dy));
 91
 92    // Copy result back to host
 93    HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(double) * m, hipMemcpyDeviceToHost));
 94
 95    std::cout << "hy" << std::endl;
 96    for(int i = 0; i < m; i++)
 97    {
 98        std::cout << hy[i] << " ";
 99    }
100    std::cout << std::endl;
101
102    // Clear hipSPARSE
103    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descr));
104    HIPSPARSE_CHECK(hipsparseDestroy(handle));
105
106    // Clear device memory
107    HIP_CHECK(hipFree(dbsrRowPtr));
108    HIP_CHECK(hipFree(dbsrColInd));
109    HIP_CHECK(hipFree(dbsrVal));
110    HIP_CHECK(hipFree(dx));
111    HIP_CHECK(hipFree(dy));
112
113    return 0;
114}

hipsparseXbsrxmv()#

hipsparseStatus_t hipsparseSbsrxmv(hipsparseHandle_t handle, hipsparseDirection_t dir, hipsparseOperation_t trans, int sizeOfMask, int mb, int nb, int nnzb, const float *alpha, const hipsparseMatDescr_t descr, const float *bsrVal, const int *bsrMaskPtr, const int *bsrRowPtr, const int *bsrEndPtr, const int *bsrColInd, int blockDim, const float *x, const float *beta, float *y)#
hipsparseStatus_t hipsparseDbsrxmv(hipsparseHandle_t handle, hipsparseDirection_t dir, hipsparseOperation_t trans, int sizeOfMask, int mb, int nb, int nnzb, const double *alpha, const hipsparseMatDescr_t descr, const double *bsrVal, const int *bsrMaskPtr, const int *bsrRowPtr, const int *bsrEndPtr, const int *bsrColInd, int blockDim, const double *x, const double *beta, double *y)#
hipsparseStatus_t hipsparseCbsrxmv(hipsparseHandle_t handle, hipsparseDirection_t dir, hipsparseOperation_t trans, int sizeOfMask, int mb, int nb, int nnzb, const hipComplex *alpha, const hipsparseMatDescr_t descr, const hipComplex *bsrVal, const int *bsrMaskPtr, const int *bsrRowPtr, const int *bsrEndPtr, const int *bsrColInd, int blockDim, const hipComplex *x, const hipComplex *beta, hipComplex *y)#
hipsparseStatus_t hipsparseZbsrxmv(hipsparseHandle_t handle, hipsparseDirection_t dir, hipsparseOperation_t trans, int sizeOfMask, int mb, int nb, int nnzb, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descr, const hipDoubleComplex *bsrVal, const int *bsrMaskPtr, const int *bsrRowPtr, const int *bsrEndPtr, const int *bsrColInd, int blockDim, const hipDoubleComplex *x, const hipDoubleComplex *beta, hipDoubleComplex *y)#

Sparse matrix vector multiplication with mask operation using BSR storage format.

hipsparseXbsrxmv multiplies the scalar \(\alpha\) with a sparse \((mb \times \text{blockDim}) \times (nb \times \text{blockDim})\) 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 == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ A^T, & \text{if trans == HIPSPARSE_OPERATION_TRANSPOSE} \\ A^H, & \text{if trans == HIPSPARSE_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, bsrRowPtr and bsr_end_ptr (both of size mb), rather the usual bsrRowPtr of size mb+1.

Deprecated:

This function is deprecated when using the CUDA backend (CUDA 12.0+) and will be removed in CUDA 13.0. This deprecation does not apply to the ROCm backend.

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 == HIPSPARSE_OPERATION_NON_TRANSPOSE is supported. Currently, blockDim == 1 is not supported.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • dir[in] matrix storage of BSR blocks.

  • trans[in] matrix operation type.

  • sizeOfMask[in] number of updated block rows of the array y. Must be non-negative and not greater than mb.

  • mb[in] number of block rows of the sparse BSR matrix. Must be non-negative.

  • nb[in] number of block columns of the sparse BSR matrix. Must be non-negative.

  • nnzb[in] number of non-zero blocks of the sparse BSR matrix. Must be non-negative.

  • alpha[in] scalar \(\alpha\).

  • descr[in] descriptor of the sparse BSR matrix. Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • bsrVal[in] array of nnzb blocks of the sparse BSR matrix.

  • bsrMaskPtr[in] array of sizeOfMask elements that give the indices of the updated block rows.

  • bsrRowPtr[in] array of mb elements that point to the start of every block row of the sparse BSR matrix.

  • bsrEndPtr[in] array of mb elements that point to the end of every block row of the sparse BSR matrix.

  • bsrColInd[in] array of nnzb elements containing the block column indices of the sparse BSR matrix.

  • blockDim[in] block dimension of the sparse BSR matrix. Must be greater than 1.

  • x[in] array of nb*blockDim elements ( \(op(A) = A\)) or mb*blockDim elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).

  • beta[in] scalar \(\beta\).

  • y[inout] array of mb*blockDim elements ( \(op(A) = A\)) or nb*blockDim elements ( \(op(A) = A^T\) or \(op(A) = A^H\)).

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, descr, alpha or beta is nullptr, mb, nb, nnzb or sizeOfMask is negative, sizeOfMask is greater than mb, blockDim is less than or equal to 1, or bsrVal, bsrMaskPtr, bsrRowPtr, bsrEndPtr, bsrColInd, x or y is nullptr.

  • HIPSPARSE_STATUS_ARCH_MISMATCH – the device is not supported.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDtrans is not HIPSPARSE_OPERATION_NON_TRANSPOSE, or hipsparseMatrixType_t is not HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXbsrsv2_zeroPivot()#

hipsparseStatus_t hipsparseXbsrsv2_zeroPivot(hipsparseHandle_t handle, bsrsv2Info_t info, int *position)#

hipsparseXbsrsv2_zeroPivot returns HIPSPARSE_STATUS_ZERO_PIVOT, if either a structural or numerical zero has been found during hipsparseXbsrsv2_analysis() or hipsparseXbsrsv2_solve() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored in position, 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 HIPSPARSE_STATUS_SUCCESS is returned instead.

Deprecated:

This function is deprecated when using the CUDA backend (CUDA 12.0+) and will be removed in CUDA 13.0. This deprecation does not apply to the ROCm backend.

Note

hipsparseXbsrsv2_zeroPivot is a blocking function. It might influence performance negatively.

Parameters:
  • handle[in] handle to the hipsparse 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:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, info or position is nullptr.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_ZERO_PIVOT – zero pivot has been found.

hipsparseXbsrsv2_bufferSize()#

hipsparseStatus_t hipsparseSbsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, float *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDbsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, double *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCbsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZbsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipDoubleComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, int *pBufferSizeInBytes)#

hipsparseXbsrsv2_bufferSize returns the size of the temporary storage buffer in bytes that is required by hipsparseXbsrsv2_analysis() and hipsparseXbsrsv2_solve(). The temporary storage buffer must be allocated by the user.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • dirA[in] matrix storage of BSR blocks.

  • transA[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.

  • descrA[in] descriptor of the sparse BSR matrix.

  • bsrSortedValA[in] array of nnzb blocks of the sparse BSR matrix.

  • bsrSortedRowPtrA[in] array of mb+1 elements that point to the start of every block row of the sparse BSR matrix.

  • bsrSortedColIndA[in] array of nnz containing the block column indices of the sparse BSR matrix.

  • blockDim[in] block dimension of the sparse BSR matrix.

  • info[out] structure that holds the information collected during the analysis step.

  • pBufferSizeInBytes[out] number of bytes of the temporary storage buffer required by hipsparseXbsrsv2_analysis() and hipsparseXbsrsv2_solve().

Return values:

hipsparseXbsrsv2_bufferSizeExt()#

hipsparseStatus_t hipsparseSbsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, float *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDbsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, double *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCbsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZbsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipDoubleComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, size_t *pBufferSizeInBytes)#

hipsparseXbsrsv2_bufferSizeExt returns the size of the temporary storage buffer in bytes that is required by hipsparseXbsrsv2_analysis() and hipsparseXbsrsv2_solve(). The temporary storage buffer must be allocated by the user.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • dirA[in] matrix storage of BSR blocks.

  • transA[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.

  • descrA[in] descriptor of the sparse BSR matrix.

  • bsrSortedValA[in] array of nnzb blocks of the sparse BSR matrix.

  • bsrSortedRowPtrA[in] array of mb+1 elements that point to the start of every block row of the sparse BSR matrix.

  • bsrSortedColIndA[in] array of nnz containing the block column indices of the sparse BSR matrix.

  • blockDim[in] block dimension of the sparse BSR matrix.

  • info[out] structure that holds the information collected during the analysis step.

  • pBufferSizeInBytes[out] number of bytes of the temporary storage buffer required by hipsparseXbsrsv2_analysis() and hipsparseXbsrsv2_solve().

Return values:

hipsparseXbsrsv2_analysis()#

hipsparseStatus_t hipsparseSbsrsv2_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, const float *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseDbsrsv2_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, const double *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseCbsrsv2_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, const hipComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseZbsrsv2_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, const hipDoubleComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#

hipsparseXbsrsv2_analysis performs the analysis step for hipsparseXbsrsv2_solve(). It is expected that this function will be executed only once for a given matrix and particular operation type.

Note

If the matrix sparsity pattern changes, the gathered information will become invalid.

Note

This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • dirA[in] matrix storage of BSR blocks.

  • transA[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.

  • descrA[in] descriptor of the sparse BSR matrix.

  • bsrSortedValA[in] array of nnzb blocks of the sparse BSR matrix.

  • bsrSortedRowPtrA[in] array of mb+1 elements that point to the start of every block row of the sparse BSR matrix.

  • bsrSortedColIndA[in] array of nnz containing the block column indices of the sparse BSR matrix.

  • blockDim[in] block dimension of the sparse BSR matrix.

  • info[out] structure that holds the information collected during the analysis step.

  • policy[in] HIPSPARSE_SOLVE_POLICY_NO_LEVEL or HIPSPARSE_SOLVE_POLICY_USE_LEVEL.

  • pBuffer[in] temporary storage buffer allocated by the user.

Return values:

hipsparseXbsrsv2_solve()#

hipsparseStatus_t hipsparseSbsrsv2_solve(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const float *alpha, const hipsparseMatDescr_t descrA, const float *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, const float *f, float *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseDbsrsv2_solve(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const double *alpha, const hipsparseMatDescr_t descrA, const double *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, const double *f, double *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseCbsrsv2_solve(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipComplex *alpha, const hipsparseMatDescr_t descrA, const hipComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, const hipComplex *f, hipComplex *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseZbsrsv2_solve(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, const hipDoubleComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, const hipDoubleComplex *f, hipDoubleComplex *x, hipsparseSolvePolicy_t policy, void *pBuffer)#

Sparse triangular solve using BSR storage format.

hipsparseXbsrsv2_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 == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ A^T, & \text{if trans == HIPSPARSE_OPERATION_TRANSPOSE} \end{array} \right. \end{split}\]

Performing the above operation requires three steps. First, the user calls hipsparseXbsrsv2_bufferSize() which will determine the size of the required temporary storage buffer. The user then allocates this buffer and calls hipsparseXbsrsv2_analysis() which will perform analysis on the sparse matrix \(op(A)\). Finally, the user completes the computation by calling hipsparseXbsrsv2_solve. The buffer size, buffer allocation, and analysis only need to be called once for a given sparse matrix \(op(A)\) while the computation stage can be repeatedly used with different \(x\) and \(y\) vectors. Once all calls to hipsparseXbsrsv2_solve are complete, the temporary buffer can be deallocated.

Solving a triangular system involves inverting the diagonal blocks. This means that if the sparse matrix is missing the diagonal block (referred to as a structural zero) or the diagonal block is not invertible (referred to as a numerical zero) then a solution is not possible. hipsparseXbsrsv2_solve tracks the location of the first zero pivot (either numerical or structural zero). The zero pivot status can be checked calling hipsparseXbsrsv2_zeroPivot(). If hipsparseXbsrsv2_zeroPivot() returns HIPSPARSE_STATUS_SUCCESS, then no zero pivot was found and therefore the matrix does not have a structural or numerical zero.

The user can specify that the sparse matrix should be interpreted as having identity blocks on the diagonal by setting the diagonal type on the descriptor descrA to HIPSPARSE_DIAG_TYPE_UNIT using hipsparseSetMatDiagType. If hipsparseDiagType_t == HIPSPARSE_DIAG_TYPE_UNIT, no zero pivot will be reported, even if the diagonal block \(A_{j,j}\) for some \(j\) is not invertible.

The sparse CSR matrix passed to hipsparseXbsrsv2_solve does not actually have to be a triangular matrix. Instead the triangular upper or lower part of the sparse matrix is solved based on hipsparseFillMode_t set on the descriptor descrA. If the fill mode is set to HIPSPARSE_FILL_MODE_LOWER, then the lower triangular matrix is solved. If the fill mode is set to HIPSPARSE_FILL_MODE_UPPER then the upper triangular matrix is solved.

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 transA == HIPSPARSE_OPERATION_NON_TRANSPOSE and transA == HIPSPARSE_OPERATION_TRANSPOSE is supported.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • dirA[in] matrix storage of BSR blocks.

  • transA[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\).

  • descrA[in] descriptor of the sparse BSR matrix.

  • bsrSortedValA[in] array of nnzb blocks of the sparse BSR matrix.

  • bsrSortedRowPtrA[in] array of mb+1 elements that point to the start of every block row of the sparse BSR matrix.

  • bsrSortedColIndA[in] array of nnz containing the block column indices of the sparse BSR matrix.

  • blockDim[in] block dimension of the sparse BSR matrix.

  • info[in] structure that holds the information collected during the analysis step.

  • f[in] array of m elements, holding the right-hand side.

  • x[out] array of m elements, holding the solution.

  • policy[in] HIPSPARSE_SOLVE_POLICY_NO_LEVEL or HIPSPARSE_SOLVE_POLICY_USE_LEVEL.

  • pBuffer[in] temporary storage buffer allocated by the user.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, mb, nnzb, blockDim, descrA, alpha, bsrSortedValA, bsrSortedRowPtrA, bsrSortedColIndA, f or x is invalid.

  • HIPSPARSE_STATUS_ARCH_MISMATCH – the device is not supported.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDtransA == HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE or hipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

  1int main(int argc, char* argv[])
  2{
  3    // hipSPARSE handle
  4    hipsparseHandle_t handle;
  5    HIPSPARSE_CHECK(hipsparseCreate(&handle));
  6
  7    // A = ( 1.0  0.0  0.0  0.0 )
  8    //     ( 2.0  3.0  0.0  0.0 )
  9    //     ( 4.0  5.0  6.0  0.0 )
 10    //     ( 7.0  0.0  8.0  9.0 )
 11    //
 12    // with bsr_dim = 2
 13    //
 14    //      -------------------
 15    //   = | 1.0 0.0 | 0.0 0.0 |
 16    //     | 2.0 3.0 | 0.0 0.0 |
 17    //      -------------------
 18    //     | 4.0 5.0 | 6.0 0.0 |
 19    //     | 7.0 0.0 | 8.0 9.0 |
 20    //      -------------------
 21
 22    // Number of rows (matrix must be square)
 23    const int m = 4;
 24
 25    // Number of block rows and block columns
 26    const int mb = 2;
 27    const int nb = 2;
 28
 29    // BSR block dimension
 30    const int bsr_dim = 2;
 31
 32    // Number of non-zero blocks
 33    const int nnzb = 3;
 34
 35    // BSR row pointers
 36    std::vector<int> hbsrRowPtr = {0, 1, 3};
 37
 38    // BSR column indices
 39    std::vector<int> hbsrColInd = {0, 0, 1};
 40
 41    // BSR values
 42    std::vector<double> hbsrVal = {1.0, 2.0, 0.0, 3.0, 4.0, 7.0, 5.0, 0.0, 6.0, 8.0, 0.0, 9.0};
 43
 44    // Storage scheme of the BSR blocks
 45    hipsparseDirection_t dir = HIPSPARSE_DIRECTION_COLUMN;
 46
 47    // Transposition of the matrix and rhs matrix
 48    hipsparseOperation_t trans = HIPSPARSE_OPERATION_NON_TRANSPOSE;
 49
 50    // Solve policy
 51    hipsparseSolvePolicy_t solve_policy = HIPSPARSE_SOLVE_POLICY_USE_LEVEL;
 52
 53    // Scalar alpha and beta
 54    double alpha = 3.7;
 55
 56    std::vector<double> hx = {1, 2, 3, 4};
 57    std::vector<double> hy(m);
 58
 59    // Offload data to device
 60    int*    dbsrRowPtr;
 61    int*    dbsrColInd;
 62    double* dbsrVal;
 63    double* dx;
 64    double* dy;
 65
 66    HIP_CHECK(hipMalloc((void**)&dbsrRowPtr, sizeof(int) * (mb + 1)));
 67    HIP_CHECK(hipMalloc((void**)&dbsrColInd, sizeof(int) * nnzb));
 68    HIP_CHECK(hipMalloc((void**)&dbsrVal, sizeof(double) * nnzb * bsr_dim * bsr_dim));
 69    HIP_CHECK(hipMalloc((void**)&dx, sizeof(double) * m));
 70    HIP_CHECK(hipMalloc((void**)&dy, sizeof(double) * m));
 71
 72    HIP_CHECK(
 73        hipMemcpy(dbsrRowPtr, hbsrRowPtr.data(), sizeof(int) * (mb + 1), hipMemcpyHostToDevice));
 74    HIP_CHECK(hipMemcpy(dbsrColInd, hbsrColInd.data(), sizeof(int) * nnzb, hipMemcpyHostToDevice));
 75    HIP_CHECK(hipMemcpy(
 76        dbsrVal, hbsrVal.data(), sizeof(double) * nnzb * bsr_dim * bsr_dim, hipMemcpyHostToDevice));
 77    HIP_CHECK(hipMemcpy(dx, hx.data(), sizeof(double) * m, hipMemcpyHostToDevice));
 78
 79    // Matrix descriptor
 80    hipsparseMatDescr_t descr;
 81    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descr));
 82
 83    // Matrix fill mode
 84    HIPSPARSE_CHECK(hipsparseSetMatFillMode(descr, HIPSPARSE_FILL_MODE_LOWER));
 85
 86    // Matrix diagonal type
 87    HIPSPARSE_CHECK(hipsparseSetMatDiagType(descr, HIPSPARSE_DIAG_TYPE_UNIT));
 88
 89    // Matrix info structure
 90    bsrsv2Info_t info;
 91    HIPSPARSE_CHECK(hipsparseCreateBsrsv2Info(&info));
 92
 93    // Obtain required buffer size
 94    int buffer_size;
 95    HIPSPARSE_CHECK(hipsparseDbsrsv2_bufferSize(handle,
 96                                                dir,
 97                                                trans,
 98                                                mb,
 99                                                nnzb,
100                                                descr,
101                                                dbsrVal,
102                                                dbsrRowPtr,
103                                                dbsrColInd,
104                                                bsr_dim,
105                                                info,
106                                                &buffer_size));
107
108    // Allocate temporary buffer
109    void* dbuffer;
110    HIP_CHECK(hipMalloc(&dbuffer, buffer_size));
111
112    // Perform analysis step
113    HIPSPARSE_CHECK(hipsparseDbsrsv2_analysis(handle,
114                                              dir,
115                                              trans,
116                                              mb,
117                                              nnzb,
118                                              descr,
119                                              dbsrVal,
120                                              dbsrRowPtr,
121                                              dbsrColInd,
122                                              bsr_dim,
123                                              info,
124                                              solve_policy,
125                                              dbuffer));
126
127    // Call dbsrsm to perform lower triangular solve LX = B
128    HIPSPARSE_CHECK(hipsparseDbsrsv2_solve(handle,
129                                           dir,
130                                           trans,
131                                           mb,
132                                           nnzb,
133                                           &alpha,
134                                           descr,
135                                           dbsrVal,
136                                           dbsrRowPtr,
137                                           dbsrColInd,
138                                           bsr_dim,
139                                           info,
140                                           dx,
141                                           dy,
142                                           solve_policy,
143                                           dbuffer));
144
145    // Check for zero pivots
146    int               pivot;
147    hipsparseStatus_t status = hipsparseXbsrsv2_zeroPivot(handle, info, &pivot);
148
149    if(status == HIPSPARSE_STATUS_ZERO_PIVOT)
150    {
151        std::cout << "Found zero pivot in matrix row " << pivot << std::endl;
152    }
153
154    // Copy results back to the host
155    HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(double) * m, hipMemcpyDeviceToHost));
156
157    std::cout << "hy" << std::endl;
158    for(int i = 0; i < m; i++)
159    {
160        std::cout << hy[i] << " ";
161    }
162    std::cout << std::endl;
163
164    // Clear hipSPARSE
165    HIPSPARSE_CHECK(hipsparseDestroyBsrsv2Info(info));
166    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descr));
167    HIPSPARSE_CHECK(hipsparseDestroy(handle));
168
169    // Clear device memory
170    HIP_CHECK(hipFree(dbsrRowPtr));
171    HIP_CHECK(hipFree(dbsrColInd));
172    HIP_CHECK(hipFree(dbsrVal));
173    HIP_CHECK(hipFree(dx));
174    HIP_CHECK(hipFree(dy));
175    HIP_CHECK(hipFree(dbuffer));
176
177    return 0;
178}

hipsparseXgemvi_bufferSize()#

hipsparseStatus_t hipsparseSgemvi_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDgemvi_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCgemvi_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZgemvi_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, int *pBufferSizeInBytes)#

hipsparseXgemvi_bufferSize returns the size of the temporary storage buffer in bytes required by hipsparseXgemvi(). The temporary storage buffer must be allocated by the user.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • transA[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.

  • pBufferSizeInBytes[out] temporary storage buffer size.

Return values:

hipsparseXgemvi()#

hipsparseStatus_t hipsparseSgemvi(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, const float *alpha, const float *A, int lda, int nnz, const float *x, const int *xInd, const float *beta, float *y, hipsparseIndexBase_t idxBase, void *pBuffer)#
hipsparseStatus_t hipsparseDgemvi(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, const double *alpha, const double *A, int lda, int nnz, const double *x, const int *xInd, const double *beta, double *y, hipsparseIndexBase_t idxBase, void *pBuffer)#
hipsparseStatus_t hipsparseCgemvi(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, const hipComplex *alpha, const hipComplex *A, int lda, int nnz, const hipComplex *x, const int *xInd, const hipComplex *beta, hipComplex *y, hipsparseIndexBase_t idxBase, void *pBuffer)#
hipsparseStatus_t hipsparseZgemvi(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, const hipDoubleComplex *alpha, const hipDoubleComplex *A, int lda, int nnz, const hipDoubleComplex *x, const int *xInd, const hipDoubleComplex *beta, hipDoubleComplex *y, hipsparseIndexBase_t idxBase, void *pBuffer)#

Dense matrix sparse vector multiplication.

hipsparseXgemvi 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
\[ op(A) = \left\{ \begin{array}{ll} A, & \text{if transA == HIPSPARSE_OPERATION_NON_TRANSPOSE} \end{array} \right. \]

Performing the above operation involves two steps. First, the user calls hipsparseXgemvi_bufferSize() in order to determine the size of the temporary storage buffer. Next, the user allocates this temporary buffer and passes it to hipsparseXgemvi to complete the computation. Once all calls to hipsparseXgemvi are complete the temporary storage buffer can be freed.

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 transA == HIPSPARSE_OPERATION_NON_TRANSPOSE is supported.

Parameters:
  • handle[in] handle to the hipsparse library context queue.

  • transA[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[in] array of nnz elements containing the values of the sparse vector

  • xInd[in] array of nnz elements containing the indices of the sparse vector

  • beta[in] scalar \(\beta\).

  • y[inout] array of m elements ( \(op(A) == A\)) or n elements ( \(op(A) == A^T\) or \(op(A) == A^H\)).

  • idxBase[in] HIPSPARSE_INDEX_BASE_ZERO or HIPSPARSE_INDEX_BASE_ONE.

  • pBuffer[in] temporary storage buffer

Return values:
 1int main(int argc, char* argv[])
 2{
 3    hipsparseOperation_t opA     = HIPSPARSE_OPERATION_NON_TRANSPOSE;
 4    hipsparseIndexBase_t idxBase = HIPSPARSE_INDEX_BASE_ZERO;
 5
 6    // Scalar alpha and beta
 7    float alpha = 1.0f;
 8    float beta  = 1.0f;
 9
10    const int m   = 4; // Number of rows of A
11    const int n   = 4; // Number of columns of A
12    const int lda = m; // leading dimension of A
13
14    // A = 1 2 3 4
15    //     5 6 7 8
16    //     2 4 6 8
17    //     4 3 2 1
18    std::vector<float> hA = {1.0f,
19                             5.0f,
20                             2.0f,
21                             4.0f,
22                             2.0f,
23                             6.0f,
24                             4.0f,
25                             3.0f,
26                             3.0f,
27                             7.0f,
28                             6.0f,
29                             2.0f,
30                             4.0f,
31                             8.0f,
32                             8.0f,
33                             1.0f};
34
35    // Sparse vector x
36    int                nnz   = 2;
37    std::vector<int>   hxInd = {0, 2};
38    std::vector<float> hx    = {10.0f, 11.0f};
39
40    // Dense vector y
41    std::vector<float> hy = {1.0f, 2.0f, 3.0f, 4.0f};
42
43    // Device data
44    float* dA = nullptr;
45    HIP_CHECK(hipMalloc((void**)&dA, sizeof(float) * m * n));
46
47    int*   dxInd = nullptr;
48    float* dx    = nullptr;
49    HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
50    HIP_CHECK(hipMalloc((void**)&dx, sizeof(float) * nnz));
51
52    float* dy = nullptr;
53    HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * m));
54
55    // Copy data from host to device
56    HIP_CHECK(hipMemcpy(dA, hA.data(), sizeof(float) * m * n, hipMemcpyHostToDevice));
57    HIP_CHECK(hipMemcpy(dxInd, hxInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
58    HIP_CHECK(hipMemcpy(dx, hx.data(), sizeof(float) * nnz, hipMemcpyHostToDevice));
59    HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(float) * m, hipMemcpyHostToDevice));
60
61    // hipSPARSE handle
62    hipsparseHandle_t handle;
63    HIPSPARSE_CHECK(hipsparseCreate(&handle));
64
65    // Call hipsparseSgemvi to perform y = alpha * A * x + beta * y
66    int bufferSize = 0;
67    HIPSPARSE_CHECK(hipsparseSgemvi_bufferSize(handle, opA, m, n, nnz, &bufferSize));
68
69    void* dbuffer = nullptr;
70    HIP_CHECK(hipMalloc((void**)&dbuffer, bufferSize));
71
72    HIPSPARSE_CHECK(hipsparseSgemvi(
73        handle, opA, m, n, &alpha, dA, lda, nnz, dx, dxInd, &beta, dy, idxBase, dbuffer));
74
75    // Copy result back to host
76    HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(float) * m, hipMemcpyDeviceToHost));
77
78    // Print the result (optional)
79    std::cout << "hy" << std::endl;
80    for(int i = 0; i < m; i++)
81    {
82        std::cout << hy[i] << " ";
83    }
84    std::cout << "" << std::endl;
85
86    // Clear hipSPARSE
87    HIPSPARSE_CHECK(hipsparseDestroy(handle));
88
89    // Clear device memory
90    HIP_CHECK(hipFree(dA));
91    HIP_CHECK(hipFree(dxInd));
92    HIP_CHECK(hipFree(dx));
93    HIP_CHECK(hipFree(dy));
94    HIP_CHECK(hipFree(dbuffer));
95
96    return 0;
97}