Sparse extra functions#

This module contains all sparse extra routines.

The sparse extra routines describe operations that manipulate sparse matrices.

hipsparseXcsrgeamNnz()#

hipsparseStatus_t hipsparseXcsrgeamNnz(hipsparseHandle_t handle, int m, int n, const hipsparseMatDescr_t descrA, int nnzA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const int *csrRowPtrB, const int *csrColIndB, const hipsparseMatDescr_t descrC, int *csrRowPtrC, int *nnzTotalDevHostPtr)#

hipsparseXcsrgeamNnz computes the total CSR non-zero elements and the CSR row offsets, that point to the start of every row of the sparse CSR matrix, of the resulting matrix \(C\). It is assumed that csrRowPtrC has been allocated with size m+1. The desired index base in the output CSR matrix is set in the hipsparseMatDescr_t. See hipsparseSetMatIndexBase().

For full code example, see hipsparseScsrgeam().

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

As indicated, nnzTotalDevHostPtr can point either to host or device memory. This is controlled by setting the pointer mode. See hipsparseSetPointerMode().

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 HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

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

  • m[in] number of rows of the sparse CSR matrix \(A\), \(B\) and \(C\). Must be non-negative.

  • n[in] number of columns of the sparse CSR matrix \(A\), \(B\) and \(C\). Must be non-negative.

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

  • nnzA[in] number of non-zero entries of the sparse CSR matrix \(A\). Must be non-negative.

  • csrRowPtrA[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(A\).

  • csrColIndA[in] array of nnzA elements containing the column indices of the sparse CSR matrix \(A\).

  • descrB[in] descriptor of the sparse CSR matrix \(B\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzB[in] number of non-zero entries of the sparse CSR matrix \(B\). Must be non-negative.

  • csrRowPtrB[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(B\).

  • csrColIndB[in] array of nnzB elements containing the column indices of the sparse CSR matrix \(B\).

  • descrC[in] descriptor of the sparse CSR matrix \(C\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • csrRowPtrC[out] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(C\).

  • nnzTotalDevHostPtr[out] pointer to the number of non-zero entries of the sparse CSR matrix \(C\). nnzTotalDevHostPtr can be a host or device pointer.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, descrA, descrB or descrC is nullptr, m, n, nnzA or nnzB is negative, or csrRowPtrA, csrColIndA, csrRowPtrB, csrColIndB, csrRowPtrC or nnzTotalDevHostPtr is nullptr.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t is not HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXcsrgeam()#

hipsparseStatus_t hipsparseScsrgeam(hipsparseHandle_t handle, int m, int n, const float *alpha, const hipsparseMatDescr_t descrA, int nnzA, const float *csrValA, const int *csrRowPtrA, const int *csrColIndA, const float *beta, const hipsparseMatDescr_t descrB, int nnzB, const float *csrValB, const int *csrRowPtrB, const int *csrColIndB, const hipsparseMatDescr_t descrC, float *csrValC, int *csrRowPtrC, int *csrColIndC)#
hipsparseStatus_t hipsparseDcsrgeam(hipsparseHandle_t handle, int m, int n, const double *alpha, const hipsparseMatDescr_t descrA, int nnzA, const double *csrValA, const int *csrRowPtrA, const int *csrColIndA, const double *beta, const hipsparseMatDescr_t descrB, int nnzB, const double *csrValB, const int *csrRowPtrB, const int *csrColIndB, const hipsparseMatDescr_t descrC, double *csrValC, int *csrRowPtrC, int *csrColIndC)#
hipsparseStatus_t hipsparseCcsrgeam(hipsparseHandle_t handle, int m, int n, const hipComplex *alpha, const hipsparseMatDescr_t descrA, int nnzA, const hipComplex *csrValA, const int *csrRowPtrA, const int *csrColIndA, const hipComplex *beta, const hipsparseMatDescr_t descrB, int nnzB, const hipComplex *csrValB, const int *csrRowPtrB, const int *csrColIndB, const hipsparseMatDescr_t descrC, hipComplex *csrValC, int *csrRowPtrC, int *csrColIndC)#
hipsparseStatus_t hipsparseZcsrgeam(hipsparseHandle_t handle, int m, int n, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, int nnzA, const hipDoubleComplex *csrValA, const int *csrRowPtrA, const int *csrColIndA, const hipDoubleComplex *beta, const hipsparseMatDescr_t descrB, int nnzB, const hipDoubleComplex *csrValB, const int *csrRowPtrB, const int *csrColIndB, const hipsparseMatDescr_t descrC, hipDoubleComplex *csrValC, int *csrRowPtrC, int *csrColIndC)#

Sparse matrix sparse matrix addition using CSR storage format.

hipsparseXcsrgeam multiplies the scalar \(\alpha\) with the sparse \(m \times n\) matrix \(A\), defined in CSR storage format, multiplies the scalar \(\beta\) with the sparse \(m \times n\) matrix \(B\), defined in CSR storage format, and adds both resulting matrices to obtain the sparse \(m \times n\) matrix \(C\), defined in CSR storage format, such that

\[ C := \alpha \cdot A + \beta \cdot B. \]

This computation involves a multi step process. First the user must allocate csrRowPtrC to have size m+1. The user then calls hipsparseXcsrgeamNnz which fills in the csrRowPtrC array as well as computes the total number of nonzeros in \(C\), nnzC. The user then allocates both arrays csrColIndC and csrValC to have size nnzC and calls hipsparseXcsrgeam to complete the computation. The desired index base in the output CSR matrix \(C\) is set in the hipsparseMatDescr_t descrC. See hipsparseSetMatIndexBase().

Note

Both scalars \(\alpha\) and \(beta\) have to be valid.

Note

Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

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.

  • m[in] number of rows of the sparse CSR matrix \(A\), \(B\) and \(C\).

  • n[in] number of columns of the sparse CSR matrix \(A\), \(B\) and \(C\).

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

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

  • nnzA[in] number of non-zero entries of the sparse CSR matrix \(A\).

  • csrValA[in] array of nnzA elements of the sparse CSR matrix \(A\).

  • csrRowPtrA[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(A\).

  • csrColIndA[in] array of nnzA elements containing the column indices of the sparse CSR matrix \(A\).

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

  • descrB[in] descriptor of the sparse CSR matrix \(B\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzB[in] number of non-zero entries of the sparse CSR matrix \(B\).

  • csrValB[in] array of nnzB elements of the sparse CSR matrix \(B\).

  • csrRowPtrB[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(B\).

  • csrColIndB[in] array of nnzB elements containing the column indices of the sparse CSR matrix \(B\).

  • descrC[in] descriptor of the sparse CSR matrix \(C\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • csrValC[out] array of elements of the sparse CSR matrix \(C\).

  • csrRowPtrC[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(C\).

  • csrColIndC[out] array of elements containing the column indices of the sparse CSR matrix \(C\).

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, n, nnzA, nnzB, alpha, descrA, csrValA, csrRowPtrA, csrColIndA, beta, descrB, csrValB, csrRowPtrB, csrColIndB, descrC, csrValC, csrRowPtrC or csrColIndC is invalid.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

  1int main(int argc, char* argv[])
  2{
  3    const int m    = 4;
  4    const int n    = 4;
  5    const int nnzA = 9;
  6    const int nnzB = 6;
  7
  8    float alpha{1.0f};
  9    float beta{1.0f};
 10
 11    // A, B, and C are m×n
 12
 13    // A
 14    // 1 0 0 2
 15    // 3 4 0 0
 16    // 5 6 7 8
 17    // 0 0 9 0
 18    std::vector<int>   hcsrRowPtrA = {0, 2, 4, 8, 9};
 19    std::vector<int>   hcsrColIndA = {0, 3, 0, 1, 0, 1, 2, 3, 2};
 20    std::vector<float> hcsrValA    = {1.0f, 2.0f, 3.0f, 4.0f, 5.0f, 6.0f, 7.0f, 8.0f, 9.0f};
 21
 22    // B
 23    // 0 1 0 0
 24    // 1 0 1 0
 25    // 0 1 0 1
 26    // 0 0 1 0
 27    std::vector<int>   hcsrRowPtrB = {0, 1, 3, 5, 6};
 28    std::vector<int>   hcsrColIndB = {1, 0, 2, 1, 3, 2};
 29    std::vector<float> hcsrValB    = {1.0f, 1.0f, 1.0f, 1.0f, 1.0f, 1.0f};
 30
 31    // Device memory management: Allocate and copy A, B
 32    int*   dcsrRowPtrA;
 33    int*   dcsrColIndA;
 34    float* dcsrValA;
 35    int*   dcsrRowPtrB;
 36    int*   dcsrColIndB;
 37    float* dcsrValB;
 38    int*   dcsrRowPtrC;
 39    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrA, (m + 1) * sizeof(int)));
 40    HIP_CHECK(hipMalloc((void**)&dcsrColIndA, nnzA * sizeof(int)));
 41    HIP_CHECK(hipMalloc((void**)&dcsrValA, nnzA * sizeof(float)));
 42    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrB, (m + 1) * sizeof(int)));
 43    HIP_CHECK(hipMalloc((void**)&dcsrColIndB, nnzB * sizeof(int)));
 44    HIP_CHECK(hipMalloc((void**)&dcsrValB, nnzB * sizeof(float)));
 45    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrC, (m + 1) * sizeof(int)));
 46
 47    HIP_CHECK(
 48        hipMemcpy(dcsrRowPtrA, hcsrRowPtrA.data(), (m + 1) * sizeof(int), hipMemcpyHostToDevice));
 49    HIP_CHECK(
 50        hipMemcpy(dcsrColIndA, hcsrColIndA.data(), nnzA * sizeof(int), hipMemcpyHostToDevice));
 51    HIP_CHECK(hipMemcpy(dcsrValA, hcsrValA.data(), nnzA * sizeof(float), hipMemcpyHostToDevice));
 52    HIP_CHECK(
 53        hipMemcpy(dcsrRowPtrB, hcsrRowPtrB.data(), (m + 1) * sizeof(int), hipMemcpyHostToDevice));
 54    HIP_CHECK(
 55        hipMemcpy(dcsrColIndB, hcsrColIndB.data(), nnzB * sizeof(int), hipMemcpyHostToDevice));
 56    HIP_CHECK(hipMemcpy(dcsrValB, hcsrValB.data(), nnzB * sizeof(float), hipMemcpyHostToDevice));
 57
 58    hipsparseHandle_t handle;
 59    HIPSPARSE_CHECK(hipsparseCreate(&handle));
 60
 61    hipsparseMatDescr_t descrA;
 62    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrA));
 63
 64    hipsparseMatDescr_t descrB;
 65    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrB));
 66
 67    hipsparseMatDescr_t descrC;
 68    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrC));
 69
 70    int nnzC;
 71    HIPSPARSE_CHECK(hipsparseXcsrgeamNnz(handle,
 72                                         m,
 73                                         n,
 74                                         descrA,
 75                                         nnzA,
 76                                         dcsrRowPtrA,
 77                                         dcsrColIndA,
 78                                         descrB,
 79                                         nnzB,
 80                                         dcsrRowPtrB,
 81                                         dcsrColIndB,
 82                                         descrC,
 83                                         dcsrRowPtrC,
 84                                         &nnzC));
 85
 86    int*   dcsrColIndC = nullptr;
 87    float* dcsrValC    = nullptr;
 88    HIP_CHECK(hipMalloc((void**)&dcsrColIndC, sizeof(int) * nnzC));
 89    HIP_CHECK(hipMalloc((void**)&dcsrValC, sizeof(float) * nnzC));
 90
 91    HIPSPARSE_CHECK(hipsparseScsrgeam(handle,
 92                                      m,
 93                                      n,
 94                                      &alpha,
 95                                      descrA,
 96                                      nnzA,
 97                                      dcsrValA,
 98                                      dcsrRowPtrA,
 99                                      dcsrColIndA,
100                                      &beta,
101                                      descrB,
102                                      nnzB,
103                                      dcsrValB,
104                                      dcsrRowPtrB,
105                                      dcsrColIndB,
106                                      descrC,
107                                      dcsrValC,
108                                      dcsrRowPtrC,
109                                      dcsrColIndC));
110
111    std::vector<int>   hcsrRowPtrC(m + 1);
112    std::vector<int>   hcsrColIndC(nnzC);
113    std::vector<float> hcsrValC(nnzC);
114
115    // Copy back to the host
116    HIP_CHECK(
117        hipMemcpy(hcsrRowPtrC.data(), dcsrRowPtrC, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
118    HIP_CHECK(
119        hipMemcpy(hcsrColIndC.data(), dcsrColIndC, sizeof(int) * nnzC, hipMemcpyDeviceToHost));
120    HIP_CHECK(hipMemcpy(hcsrValC.data(), dcsrValC, sizeof(float) * nnzC, hipMemcpyDeviceToHost));
121
122    std::cout << "C" << std::endl;
123    for(int i = 0; i < m; i++)
124    {
125        int start = hcsrRowPtrC[i];
126        int end   = hcsrRowPtrC[i + 1];
127
128        std::vector<float> temp(n, 0.0f);
129        for(int j = start; j < end; j++)
130        {
131            temp[hcsrColIndC[j]] = hcsrValC[j];
132        }
133
134        for(int j = 0; j < n; j++)
135        {
136            std::cout << temp[j] << " ";
137        }
138        std::cout << std::endl;
139    }
140    std::cout << std::endl;
141
142    HIP_CHECK(hipFree(dcsrRowPtrA));
143    HIP_CHECK(hipFree(dcsrColIndA));
144    HIP_CHECK(hipFree(dcsrValA));
145    HIP_CHECK(hipFree(dcsrRowPtrB));
146    HIP_CHECK(hipFree(dcsrColIndB));
147    HIP_CHECK(hipFree(dcsrValB));
148    HIP_CHECK(hipFree(dcsrRowPtrC));
149    HIP_CHECK(hipFree(dcsrColIndC));
150    HIP_CHECK(hipFree(dcsrValC));
151
152    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrA));
153    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrB));
154    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrC));
155    HIPSPARSE_CHECK(hipsparseDestroy(handle));
156
157    return 0;
158}

hipsparseXcsrgeam2_bufferSizeExt()#

hipsparseStatus_t hipsparseScsrgeam2_bufferSizeExt(hipsparseHandle_t handle, int m, int n, const float *alpha, const hipsparseMatDescr_t descrA, int nnzA, const float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const float *beta, const hipsparseMatDescr_t descrB, int nnzB, const float *csrSortedValB, const int *csrSortedRowPtrB, const int *csrSortedColIndB, const hipsparseMatDescr_t descrC, const float *csrSortedValC, const int *csrSortedRowPtrC, const int *csrSortedColIndC, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDcsrgeam2_bufferSizeExt(hipsparseHandle_t handle, int m, int n, const double *alpha, const hipsparseMatDescr_t descrA, int nnzA, const double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const double *beta, const hipsparseMatDescr_t descrB, int nnzB, const double *csrSortedValB, const int *csrSortedRowPtrB, const int *csrSortedColIndB, const hipsparseMatDescr_t descrC, const double *csrSortedValC, const int *csrSortedRowPtrC, const int *csrSortedColIndC, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCcsrgeam2_bufferSizeExt(hipsparseHandle_t handle, int m, int n, const hipComplex *alpha, const hipsparseMatDescr_t descrA, int nnzA, const hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const hipComplex *beta, const hipsparseMatDescr_t descrB, int nnzB, const hipComplex *csrSortedValB, const int *csrSortedRowPtrB, const int *csrSortedColIndB, const hipsparseMatDescr_t descrC, const hipComplex *csrSortedValC, const int *csrSortedRowPtrC, const int *csrSortedColIndC, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZcsrgeam2_bufferSizeExt(hipsparseHandle_t handle, int m, int n, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, int nnzA, const hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const hipDoubleComplex *beta, const hipsparseMatDescr_t descrB, int nnzB, const hipDoubleComplex *csrSortedValB, const int *csrSortedRowPtrB, const int *csrSortedColIndB, const hipsparseMatDescr_t descrC, const hipDoubleComplex *csrSortedValC, const int *csrSortedRowPtrC, const int *csrSortedColIndC, size_t *pBufferSizeInBytes)#

hipsparseXcsrgeam2_bufferSizeExt returns the size of the temporary storage buffer in bytes that is required by hipsparseXcsrgeam2Nnz() and hipsparseXcsrgeam2(). The temporary storage buffer must be allocated by the user.

Note

Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

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

  • m[in] number of rows of the sparse CSR matrix \(A\), \(B\) and \(C\).

  • n[in] number of columns of the sparse CSR matrix \(A\), \(B\) and \(C\).

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

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

  • nnzA[in] number of non-zero entries of the sparse CSR matrix \(A\).

  • csrSortedValA[in] array of nnzA elements of the sparse CSR matrix \(A\).

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

  • csrSortedColIndA[in] array of nnzA elements containing the column indices of the sparse CSR matrix \(A\).

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

  • descrB[in] descriptor of the sparse CSR matrix \(B\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzB[in] number of non-zero entries of the sparse CSR matrix \(B\).

  • csrSortedValB[in] array of nnzB elements of the sparse CSR matrix \(B\).

  • csrSortedRowPtrB[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(B\).

  • csrSortedColIndB[in] array of nnzB elements containing the column indices of the sparse CSR matrix \(B\).

  • descrC[in] descriptor of the sparse CSR matrix \(C\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • csrSortedValC[out] array of elements of the sparse CSR matrix \(C\).

  • csrSortedRowPtrC[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(C\).

  • csrSortedColIndC[out] array of elements containing the column indices of the sparse CSR matrix \(C\).

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, n, nnzA, nnzB, alpha, descrA, csrSortedValA, csrSortedRowPtrA, csrSortedColIndA, beta, descrB, csrSortedValB, csrSortedRowPtrB, csrSortedColIndB, descrC, csrSortedValC, csrSortedRowPtrC, csrSortedColIndC, or pBufferSizeInBytes is invalid.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXcsrgeam2Nnz()#

hipsparseStatus_t hipsparseXcsrgeam2Nnz(hipsparseHandle_t handle, int m, int n, const hipsparseMatDescr_t descrA, int nnzA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const hipsparseMatDescr_t descrB, int nnzB, const int *csrSortedRowPtrB, const int *csrSortedColIndB, const hipsparseMatDescr_t descrC, int *csrSortedRowPtrC, int *nnzTotalDevHostPtr, void *workspace)#

hipsparseXcsrgeam2Nnz computes the total CSR non-zero elements and the CSR row offsets, that point to the start of every row of the sparse CSR matrix, of the resulting matrix \(C\). It is assumed that csrRowPtrC has been allocated with size m+1. The required buffer size can be obtained by hipsparseXcsrgeam2_bufferSizeExt(). The desired index base in the output CSR matrix \(C\) is set in the hipsparseMatDescr_t descrC. See hipsparseSetMatIndexBase().

Note

As indicated, nnzTotalDevHostPtr can point either to host or device memory. This is controlled by setting the pointer mode. See hipsparseSetPointerMode().

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 HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

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

  • m[in] number of rows of the sparse CSR matrix \(A\), \(B\) and \(C\).

  • n[in] number of columns of the sparse CSR matrix \(A\), \(B\) and \(C\).

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

  • nnzA[in] number of non-zero entries of the sparse CSR matrix \(A\).

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

  • csrSortedColIndA[in] array of nnzA elements containing the column indices of the sparse CSR matrix \(A\).

  • descrB[in] descriptor of the sparse CSR matrix \(B\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzB[in] number of non-zero entries of the sparse CSR matrix \(B\).

  • csrSortedRowPtrB[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(B\).

  • csrSortedColIndB[in] array of nnzB elements containing the column indices of the sparse CSR matrix \(B\).

  • descrC[in] descriptor of the sparse CSR matrix \(C\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • csrSortedRowPtrC[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(C\).

  • nnzTotalDevHostPtr[out] pointer to the number of non-zero entries of the sparse CSR matrix \(C\). nnzTotalDevHostPtr can be a host or device pointer.

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, n, nnzA, nnzB, descrA, csrSortedRowPtrA, csrSortedColIndA, descrB, csrSortedRowPtrB, csrSortedColIndB, descrC, csrSortedRowPtrC or nnzTotalDevHostPtr is invalid.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXcsrgeam2()#

hipsparseStatus_t hipsparseScsrgeam2(hipsparseHandle_t handle, int m, int n, const float *alpha, const hipsparseMatDescr_t descrA, int nnzA, const float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const float *beta, const hipsparseMatDescr_t descrB, int nnzB, const float *csrSortedValB, const int *csrSortedRowPtrB, const int *csrSortedColIndB, const hipsparseMatDescr_t descrC, float *csrSortedValC, int *csrSortedRowPtrC, int *csrSortedColIndC, void *pBuffer)#
hipsparseStatus_t hipsparseDcsrgeam2(hipsparseHandle_t handle, int m, int n, const double *alpha, const hipsparseMatDescr_t descrA, int nnzA, const double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const double *beta, const hipsparseMatDescr_t descrB, int nnzB, const double *csrSortedValB, const int *csrSortedRowPtrB, const int *csrSortedColIndB, const hipsparseMatDescr_t descrC, double *csrSortedValC, int *csrSortedRowPtrC, int *csrSortedColIndC, void *pBuffer)#
hipsparseStatus_t hipsparseCcsrgeam2(hipsparseHandle_t handle, int m, int n, const hipComplex *alpha, const hipsparseMatDescr_t descrA, int nnzA, const hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const hipComplex *beta, const hipsparseMatDescr_t descrB, int nnzB, const hipComplex *csrSortedValB, const int *csrSortedRowPtrB, const int *csrSortedColIndB, const hipsparseMatDescr_t descrC, hipComplex *csrSortedValC, int *csrSortedRowPtrC, int *csrSortedColIndC, void *pBuffer)#
hipsparseStatus_t hipsparseZcsrgeam2(hipsparseHandle_t handle, int m, int n, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, int nnzA, const hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const hipDoubleComplex *beta, const hipsparseMatDescr_t descrB, int nnzB, const hipDoubleComplex *csrSortedValB, const int *csrSortedRowPtrB, const int *csrSortedColIndB, const hipsparseMatDescr_t descrC, hipDoubleComplex *csrSortedValC, int *csrSortedRowPtrC, int *csrSortedColIndC, void *pBuffer)#

Sparse matrix sparse matrix addition using CSR storage format.

hipsparseXcsrgeam2 multiplies the scalar \(\alpha\) with the sparse \(m \times n\) matrix \(A\), defined in CSR storage format, multiplies the scalar \(\beta\) with the sparse \(m \times n\) matrix \(B\), defined in CSR storage format, and adds both resulting matrices to obtain the sparse \(m \times n\) matrix \(C\), defined in CSR storage format, such that

\[ C := \alpha \cdot A + \beta \cdot B. \]

This computation involves a multi step process. First the user must call hipsparseXcsrgeam2_bufferSizeExt() in order to determine the required user allocated temporary buffer size. The user then allocates this buffer and also allocates csrRowPtrC to have size m+1. Both the temporary storage buffer and csrRowPtrC array are then passed to hipsparseXcsrgeam2Nnz which fills in the csrRowPtrC array as well as computes the total number of nonzeros in C, nnzC. The user then allocates both arrays csrColIndC and csrValC to have size nnzC and calls hipsparseXcsrgeam2 to complete the computation. The desired index base in the output CSR matrix C is set in the hipsparseMatDescr_t descrC. See hipsparseSetMatIndexBase().

Note

Both scalars \(\alpha\) and \(beta\) have to be valid.

Note

Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

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.

  • m[in] number of rows of the sparse CSR matrix \(A\), \(B\) and \(C\).

  • n[in] number of columns of the sparse CSR matrix \(A\), \(B\) and \(C\).

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

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

  • nnzA[in] number of non-zero entries of the sparse CSR matrix \(A\).

  • csrSortedValA[in] array of nnzA elements of the sparse CSR matrix \(A\).

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

  • csrSortedColIndA[in] array of nnzA elements containing the column indices of the sparse CSR matrix \(A\).

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

  • descrB[in] descriptor of the sparse CSR matrix \(B\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzB[in] number of non-zero entries of the sparse CSR matrix \(B\).

  • csrSortedValB[in] array of nnzB elements of the sparse CSR matrix \(B\).

  • csrSortedRowPtrB[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(B\).

  • csrSortedColIndB[in] array of nnzB elements containing the column indices of the sparse CSR matrix \(B\).

  • descrC[in] descriptor of the sparse CSR matrix \(C\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • csrSortedValC[out] array of elements of the sparse CSR matrix \(C\).

  • csrSortedRowPtrC[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(C\).

  • csrSortedColIndC[out] array of elements containing the column indices of the sparse CSR matrix \(C\).

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, n, nnzA, nnzB, alpha, descrA, csrSortedValA, csrSortedRowPtrA, csrSortedColIndA, beta, descrB, csrSortedValB, csrSortedRowPtrB, csrSortedColIndB, descrC, csrSortedValC, csrSortedRowPtrC, csrSortedColIndC or pBuffer is invalid.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

  1int main(int argc, char* argv[])
  2{
  3    const int m    = 4;
  4    const int n    = 4;
  5    const int nnzA = 9;
  6    const int nnzB = 6;
  7
  8    float alpha{1.0f};
  9    float beta{1.0f};
 10
 11    // A, B, and C are m×n
 12
 13    // A
 14    // 1 0 0 2
 15    // 3 4 0 0
 16    // 5 6 7 8
 17    // 0 0 9 0
 18    std::vector<int>   hcsrRowPtrA = {0, 2, 4, 8, 9};
 19    std::vector<int>   hcsrColIndA = {0, 3, 0, 1, 0, 1, 2, 3, 2};
 20    std::vector<float> hcsrValA    = {1.0f, 2.0f, 3.0f, 4.0f, 5.0f, 6.0f, 7.0f, 8.0f, 9.0f};
 21
 22    // B
 23    // 0 1 0 0
 24    // 1 0 1 0
 25    // 0 1 0 1
 26    // 0 0 1 0
 27    std::vector<int>   hcsrRowPtrB = {0, 1, 3, 5, 6};
 28    std::vector<int>   hcsrColIndB = {1, 0, 2, 1, 3, 2};
 29    std::vector<float> hcsrValB    = {1.0f, 1.0f, 1.0f, 1.0f, 1.0f, 1.0f};
 30
 31    // Device memory management: Allocate and copy A, B
 32    int*   dcsrRowPtrA;
 33    int*   dcsrColIndA;
 34    float* dcsrValA;
 35    int*   dcsrRowPtrB;
 36    int*   dcsrColIndB;
 37    float* dcsrValB;
 38    int*   dcsrRowPtrC;
 39    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrA, (m + 1) * sizeof(int)));
 40    HIP_CHECK(hipMalloc((void**)&dcsrColIndA, nnzA * sizeof(int)));
 41    HIP_CHECK(hipMalloc((void**)&dcsrValA, nnzA * sizeof(float)));
 42    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrB, (m + 1) * sizeof(int)));
 43    HIP_CHECK(hipMalloc((void**)&dcsrColIndB, nnzB * sizeof(int)));
 44    HIP_CHECK(hipMalloc((void**)&dcsrValB, nnzB * sizeof(float)));
 45    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrC, (m + 1) * sizeof(int)));
 46
 47    HIP_CHECK(
 48        hipMemcpy(dcsrRowPtrA, hcsrRowPtrA.data(), (m + 1) * sizeof(int), hipMemcpyHostToDevice));
 49    HIP_CHECK(
 50        hipMemcpy(dcsrColIndA, hcsrColIndA.data(), nnzA * sizeof(int), hipMemcpyHostToDevice));
 51    HIP_CHECK(hipMemcpy(dcsrValA, hcsrValA.data(), nnzA * sizeof(float), hipMemcpyHostToDevice));
 52    HIP_CHECK(
 53        hipMemcpy(dcsrRowPtrB, hcsrRowPtrB.data(), (m + 1) * sizeof(int), hipMemcpyHostToDevice));
 54    HIP_CHECK(
 55        hipMemcpy(dcsrColIndB, hcsrColIndB.data(), nnzB * sizeof(int), hipMemcpyHostToDevice));
 56    HIP_CHECK(hipMemcpy(dcsrValB, hcsrValB.data(), nnzB * sizeof(float), hipMemcpyHostToDevice));
 57
 58    hipsparseHandle_t handle;
 59    HIPSPARSE_CHECK(hipsparseCreate(&handle));
 60
 61    hipsparseMatDescr_t descrA;
 62    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrA));
 63
 64    hipsparseMatDescr_t descrB;
 65    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrB));
 66
 67    hipsparseMatDescr_t descrC;
 68    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrC));
 69
 70    size_t bufferSize;
 71    HIPSPARSE_CHECK(hipsparseScsrgeam2_bufferSizeExt(handle,
 72                                                     m,
 73                                                     n,
 74                                                     &alpha,
 75                                                     descrA,
 76                                                     nnzA,
 77                                                     dcsrValA,
 78                                                     dcsrRowPtrA,
 79                                                     dcsrColIndA,
 80                                                     &beta,
 81                                                     descrB,
 82                                                     nnzB,
 83                                                     dcsrValB,
 84                                                     dcsrRowPtrB,
 85                                                     dcsrColIndB,
 86                                                     descrC,
 87                                                     nullptr,
 88                                                     dcsrRowPtrC,
 89                                                     nullptr,
 90                                                     &bufferSize));
 91
 92    void* dbuffer = nullptr;
 93    HIP_CHECK(hipMalloc((void**)&dbuffer, bufferSize));
 94
 95    int nnzC;
 96    HIPSPARSE_CHECK(hipsparseXcsrgeam2Nnz(handle,
 97                                          m,
 98                                          n,
 99                                          descrA,
100                                          nnzA,
101                                          dcsrRowPtrA,
102                                          dcsrColIndA,
103                                          descrB,
104                                          nnzB,
105                                          dcsrRowPtrB,
106                                          dcsrColIndB,
107                                          descrC,
108                                          dcsrRowPtrC,
109                                          &nnzC,
110                                          dbuffer));
111
112    int*   dcsrColIndC = nullptr;
113    float* dcsrValC    = nullptr;
114    HIP_CHECK(hipMalloc((void**)&dcsrColIndC, sizeof(int) * nnzC));
115    HIP_CHECK(hipMalloc((void**)&dcsrValC, sizeof(float) * nnzC));
116
117    HIPSPARSE_CHECK(hipsparseScsrgeam2(handle,
118                                       m,
119                                       n,
120                                       &alpha,
121                                       descrA,
122                                       nnzA,
123                                       dcsrValA,
124                                       dcsrRowPtrA,
125                                       dcsrColIndA,
126                                       &beta,
127                                       descrB,
128                                       nnzB,
129                                       dcsrValB,
130                                       dcsrRowPtrB,
131                                       dcsrColIndB,
132                                       descrC,
133                                       dcsrValC,
134                                       dcsrRowPtrC,
135                                       dcsrColIndC,
136                                       dbuffer));
137
138    std::vector<int>   hcsrRowPtrC(m + 1);
139    std::vector<int>   hcsrColIndC(nnzC);
140    std::vector<float> hcsrValC(nnzC);
141
142    // Copy back to the host
143    HIP_CHECK(
144        hipMemcpy(hcsrRowPtrC.data(), dcsrRowPtrC, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
145    HIP_CHECK(
146        hipMemcpy(hcsrColIndC.data(), dcsrColIndC, sizeof(int) * nnzC, hipMemcpyDeviceToHost));
147    HIP_CHECK(hipMemcpy(hcsrValC.data(), dcsrValC, sizeof(float) * nnzC, hipMemcpyDeviceToHost));
148
149    std::cout << "C" << std::endl;
150    for(int i = 0; i < m; i++)
151    {
152        int start = hcsrRowPtrC[i];
153        int end   = hcsrRowPtrC[i + 1];
154
155        std::vector<float> temp(n, 0.0f);
156        for(int j = start; j < end; j++)
157        {
158            temp[hcsrColIndC[j]] = hcsrValC[j];
159        }
160
161        for(int j = 0; j < n; j++)
162        {
163            std::cout << temp[j] << " ";
164        }
165        std::cout << std::endl;
166    }
167    std::cout << std::endl;
168
169    HIP_CHECK(hipFree(dcsrRowPtrA));
170    HIP_CHECK(hipFree(dcsrColIndA));
171    HIP_CHECK(hipFree(dcsrValA));
172    HIP_CHECK(hipFree(dcsrRowPtrB));
173    HIP_CHECK(hipFree(dcsrColIndB));
174    HIP_CHECK(hipFree(dcsrValB));
175    HIP_CHECK(hipFree(dcsrRowPtrC));
176    HIP_CHECK(hipFree(dcsrColIndC));
177    HIP_CHECK(hipFree(dcsrValC));
178
179    HIP_CHECK(hipFree(dbuffer));
180
181    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrA));
182    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrB));
183    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrC));
184    HIPSPARSE_CHECK(hipsparseDestroy(handle));
185
186    return 0;
187}

hipsparseXcsrgemmNnz()#

hipsparseStatus_t hipsparseXcsrgemmNnz(hipsparseHandle_t handle, hipsparseOperation_t transA, hipsparseOperation_t transB, int m, int n, int k, const hipsparseMatDescr_t descrA, int nnzA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const int *csrRowPtrB, const int *csrColIndB, const hipsparseMatDescr_t descrC, int *csrRowPtrC, int *nnzTotalDevHostPtr)#

hipsparseXcsrgemmNnz computes the total CSR non-zero elements and the CSR row offsets, that point to the start of every row of the sparse CSR matrix, of the resulting multiplied matrix \(C\). It is assumed that csrRowPtrC has been allocated with size m+1. The desired index base in the output CSR matrix \(C\) is set in the hipsparseMatDescr_t descrC. See hipsparseSetMatIndexBase().

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

As indicated, nnzTotalDevHostPtr can point either to host or device memory. This is controlled by setting the pointer mode. See hipsparseSetPointerMode().

Note

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

Note

Please note, that for matrix products with more than 8192 intermediate products per row, additional temporary storage buffer is allocated by the algorithm.

Note

Currently, only transA == transB == HIPSPARSE_OPERATION_NON_TRANSPOSE is supported.

Note

Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

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

  • transA[in] matrix \(A\) operation type.

  • transB[in] matrix \(B\) operation type.

  • m[in] number of rows of the sparse CSR matrix \(op(A)\) and \(C\). Must be non-negative.

  • n[in] number of columns of the sparse CSR matrix \(op(B)\) and \(C\). Must be non-negative.

  • k[in] number of columns of the sparse CSR matrix \(op(A)\) and number of rows of the sparse CSR matrix \(op(B)\). Must be non-negative.

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

  • nnzA[in] number of non-zero entries of the sparse CSR matrix \(A\). Must be non-negative.

  • csrRowPtrA[in] array of m+1 elements ( \(op(A) == A\), k+1 otherwise) that point to the start of every row of the sparse CSR matrix \(op(A)\).

  • csrColIndA[in] array of nnzA elements containing the column indices of the sparse CSR matrix \(A\).

  • descrB[in] descriptor of the sparse CSR matrix \(B\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzB[in] number of non-zero entries of the sparse CSR matrix \(B\). Must be non-negative.

  • csrRowPtrB[in] array of k+1 elements ( \(op(B) == B\), m+1 otherwise) that point to the start of every row of the sparse CSR matrix \(op(B)\).

  • csrColIndB[in] array of nnzB elements containing the column indices of the sparse CSR matrix \(B\).

  • descrC[in] descriptor of the sparse CSR matrix \(C\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • csrRowPtrC[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(C\).

  • nnzTotalDevHostPtr[inout] pointer to the number of non-zero entries of the sparse CSR matrix \(C\). nnzTotalDevHostPtr can be a host or device pointer.

Return values:

hipsparseXcsrgemm()#

hipsparseStatus_t hipsparseScsrgemm(hipsparseHandle_t handle, hipsparseOperation_t transA, hipsparseOperation_t transB, int m, int n, int k, const hipsparseMatDescr_t descrA, int nnzA, const float *csrValA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const float *csrValB, const int *csrRowPtrB, const int *csrColIndB, const hipsparseMatDescr_t descrC, float *csrValC, const int *csrRowPtrC, int *csrColIndC)#
hipsparseStatus_t hipsparseDcsrgemm(hipsparseHandle_t handle, hipsparseOperation_t transA, hipsparseOperation_t transB, int m, int n, int k, const hipsparseMatDescr_t descrA, int nnzA, const double *csrValA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const double *csrValB, const int *csrRowPtrB, const int *csrColIndB, const hipsparseMatDescr_t descrC, double *csrValC, const int *csrRowPtrC, int *csrColIndC)#
hipsparseStatus_t hipsparseCcsrgemm(hipsparseHandle_t handle, hipsparseOperation_t transA, hipsparseOperation_t transB, int m, int n, int k, const hipsparseMatDescr_t descrA, int nnzA, const hipComplex *csrValA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const hipComplex *csrValB, const int *csrRowPtrB, const int *csrColIndB, const hipsparseMatDescr_t descrC, hipComplex *csrValC, const int *csrRowPtrC, int *csrColIndC)#
hipsparseStatus_t hipsparseZcsrgemm(hipsparseHandle_t handle, hipsparseOperation_t transA, hipsparseOperation_t transB, int m, int n, int k, const hipsparseMatDescr_t descrA, int nnzA, const hipDoubleComplex *csrValA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const hipDoubleComplex *csrValB, const int *csrRowPtrB, const int *csrColIndB, const hipsparseMatDescr_t descrC, hipDoubleComplex *csrValC, const int *csrRowPtrC, int *csrColIndC)#

Sparse matrix sparse matrix multiplication using CSR storage format.

hipsparseXcsrgemm multiplies the sparse \(m \times k\) matrix \(op(A)\), defined in CSR storage format with the sparse \(k \times n\) matrix \(op(B)\), defined in CSR storage format, and stores the result in the sparse \(m \times n\) matrix \(C\), defined in CSR storage format, such that

\[ C := op(A) \cdot op(B), \]
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}\]
and
\[\begin{split} op(B) = \left\{ \begin{array}{ll} B, & \text{if transB == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ B^T, & \text{if transB == HIPSPARSE_OPERATION_TRANSPOSE} \\ B^H, & \text{if transB == HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE} \end{array} \right. \end{split}\]

This computation involves a multi step process. First the user must allocate csrRowPtrC to have size m+1. The user then calls hipsparseXcsrgemmNnz which fills in the csrRowPtrC array as well as computes the total number of nonzeros in C, nnzC. The user then allocates both arrays csrColIndC and csrValC to have size nnzC and calls hipsparseXcsrgemm to complete the computation. The desired index base in the output CSR matrix C is set in the hipsparseMatDescr_t descrC. See hipsparseSetMatIndexBase().

Note

Currently, only transA == HIPSPARSE_OPERATION_NON_TRANSPOSE is supported.

Note

Currently, only transB == HIPSPARSE_OPERATION_NON_TRANSPOSE is supported.

Note

Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

Note

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

Note

Please note, that for matrix products with more than 4096 non-zero entries per row, additional temporary storage buffer is allocated by the algorithm.

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

  • transA[in] matrix \(A\) operation type.

  • transB[in] matrix \(B\) operation type.

  • m[in] number of rows of the sparse CSR matrix \(op(A)\) and \(C\).

  • n[in] number of columns of the sparse CSR matrix \(op(B)\) and \(C\).

  • k[in] number of columns of the sparse CSR matrix \(op(A)\) and number of rows of the sparse CSR matrix \(op(B)\).

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

  • nnzA[in] number of non-zero entries of the sparse CSR matrix \(A\).

  • csrValA[in] array of nnzA elements of the sparse CSR matrix \(A\).

  • csrRowPtrA[in] array of m+1 elements ( \(op(A) == A\), k+1 otherwise) that point to the start of every row of the sparse CSR matrix \(op(A)\).

  • csrColIndA[in] array of nnzA elements containing the column indices of the sparse CSR matrix \(A\).

  • descrB[in] descriptor of the sparse CSR matrix \(B\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzB[in] number of non-zero entries of the sparse CSR matrix \(B\).

  • csrValB[in] array of nnzB elements of the sparse CSR matrix \(B\).

  • csrRowPtrB[in] array of k+1 elements ( \(op(B) == B\), m+1 otherwise) that point to the start of every row of the sparse CSR matrix \(op(B)\).

  • csrColIndB[in] array of nnzB elements containing the column indices of the sparse CSR matrix \(B\).

  • descrC[in] descriptor of the sparse CSR matrix \(C\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • csrValC[out] array of nnzC elements of the sparse CSR matrix \(C\).

  • csrRowPtrC[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(C\).

  • csrColIndC[out] array of nnzC elements containing the column indices of the sparse CSR matrix \(C\).

Return values:
  1int main(int argc, char* argv[])
  2{
  3    const int m    = 4;
  4    const int k    = 3;
  5    const int n    = 2;
  6    const int nnzA = 7;
  7    const int nnzB = 3;
  8
  9    hipsparseOperation_t transA = HIPSPARSE_OPERATION_NON_TRANSPOSE;
 10    hipsparseOperation_t transB = HIPSPARSE_OPERATION_NON_TRANSPOSE;
 11
 12    // A, B, and C are mxk, kxn, and m×n
 13
 14    // A
 15    // 1 0 0
 16    // 3 4 0
 17    // 5 6 7
 18    // 0 0 9
 19    std::vector<int>   hcsrRowPtrA = {0, 1, 3, 6, 7};
 20    std::vector<int>   hcsrColIndA = {0, 0, 1, 0, 1, 2, 2};
 21    std::vector<float> hcsrValA    = {1.0f, 3.0f, 4.0f, 5.0f, 6.0f, 7.0f, 9.0f};
 22
 23    // B
 24    // 0 1
 25    // 1 0
 26    // 0 1
 27    std::vector<int>   hcsrRowPtrB = {0, 1, 2, 3};
 28    std::vector<int>   hcsrColIndB = {1, 0, 1};
 29    std::vector<float> hcsrValB    = {1.0f, 1.0f, 1.0f};
 30
 31    // Device memory management: Allocate and copy A, B
 32    int*   dcsrRowPtrA;
 33    int*   dcsrColIndA;
 34    float* dcsrValA;
 35    int*   dcsrRowPtrB;
 36    int*   dcsrColIndB;
 37    float* dcsrValB;
 38    int*   dcsrRowPtrC;
 39    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrA, (m + 1) * sizeof(int)));
 40    HIP_CHECK(hipMalloc((void**)&dcsrColIndA, nnzA * sizeof(int)));
 41    HIP_CHECK(hipMalloc((void**)&dcsrValA, nnzA * sizeof(float)));
 42    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrB, (m + 1) * sizeof(int)));
 43    HIP_CHECK(hipMalloc((void**)&dcsrColIndB, nnzB * sizeof(int)));
 44    HIP_CHECK(hipMalloc((void**)&dcsrValB, nnzB * sizeof(float)));
 45    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrC, (m + 1) * sizeof(int)));
 46
 47    HIP_CHECK(
 48        hipMemcpy(dcsrRowPtrA, hcsrRowPtrA.data(), (m + 1) * sizeof(int), hipMemcpyHostToDevice));
 49    HIP_CHECK(
 50        hipMemcpy(dcsrColIndA, hcsrColIndA.data(), nnzA * sizeof(int), hipMemcpyHostToDevice));
 51    HIP_CHECK(hipMemcpy(dcsrValA, hcsrValA.data(), nnzA * sizeof(float), hipMemcpyHostToDevice));
 52    HIP_CHECK(
 53        hipMemcpy(dcsrRowPtrB, hcsrRowPtrB.data(), (m + 1) * sizeof(int), hipMemcpyHostToDevice));
 54    HIP_CHECK(
 55        hipMemcpy(dcsrColIndB, hcsrColIndB.data(), nnzB * sizeof(int), hipMemcpyHostToDevice));
 56    HIP_CHECK(hipMemcpy(dcsrValB, hcsrValB.data(), nnzB * sizeof(float), hipMemcpyHostToDevice));
 57
 58    hipsparseHandle_t handle;
 59    HIPSPARSE_CHECK(hipsparseCreate(&handle));
 60
 61    hipsparseMatDescr_t descrA;
 62    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrA));
 63
 64    hipsparseMatDescr_t descrB;
 65    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrB));
 66
 67    hipsparseMatDescr_t descrC;
 68    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrC));
 69
 70    int nnzC;
 71    HIPSPARSE_CHECK(hipsparseXcsrgemmNnz(handle,
 72                                         transA,
 73                                         transB,
 74                                         m,
 75                                         n,
 76                                         k,
 77                                         descrA,
 78                                         nnzA,
 79                                         dcsrRowPtrA,
 80                                         dcsrColIndA,
 81                                         descrB,
 82                                         nnzB,
 83                                         dcsrRowPtrB,
 84                                         dcsrColIndB,
 85                                         descrC,
 86                                         dcsrRowPtrC,
 87                                         &nnzC));
 88
 89    int*   dcsrColIndC = nullptr;
 90    float* dcsrValC    = nullptr;
 91    HIP_CHECK(hipMalloc((void**)&dcsrColIndC, sizeof(int) * nnzC));
 92    HIP_CHECK(hipMalloc((void**)&dcsrValC, sizeof(float) * nnzC));
 93
 94    HIPSPARSE_CHECK(hipsparseScsrgemm(handle,
 95                                      transA,
 96                                      transB,
 97                                      m,
 98                                      n,
 99                                      k,
100                                      descrA,
101                                      nnzA,
102                                      dcsrValA,
103                                      dcsrRowPtrA,
104                                      dcsrColIndA,
105                                      descrB,
106                                      nnzB,
107                                      dcsrValB,
108                                      dcsrRowPtrB,
109                                      dcsrColIndB,
110                                      descrC,
111                                      dcsrValC,
112                                      dcsrRowPtrC,
113                                      dcsrColIndC));
114
115    std::vector<int>   hcsrRowPtrC(m + 1);
116    std::vector<int>   hcsrColIndC(nnzC);
117    std::vector<float> hcsrValC(nnzC);
118
119    // Copy back to the host
120    HIP_CHECK(
121        hipMemcpy(hcsrRowPtrC.data(), dcsrRowPtrC, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
122    HIP_CHECK(
123        hipMemcpy(hcsrColIndC.data(), dcsrColIndC, sizeof(int) * nnzC, hipMemcpyDeviceToHost));
124    HIP_CHECK(hipMemcpy(hcsrValC.data(), dcsrValC, sizeof(float) * nnzC, hipMemcpyDeviceToHost));
125
126    std::cout << "C" << std::endl;
127    for(int i = 0; i < m; i++)
128    {
129        int start = hcsrRowPtrC[i];
130        int end   = hcsrRowPtrC[i + 1];
131
132        std::vector<float> temp(n, 0.0f);
133        for(int j = start; j < end; j++)
134        {
135            temp[hcsrColIndC[j]] = hcsrValC[j];
136        }
137
138        for(int j = 0; j < n; j++)
139        {
140            std::cout << temp[j] << " ";
141        }
142        std::cout << std::endl;
143    }
144    std::cout << std::endl;
145
146    HIP_CHECK(hipFree(dcsrRowPtrA));
147    HIP_CHECK(hipFree(dcsrColIndA));
148    HIP_CHECK(hipFree(dcsrValA));
149    HIP_CHECK(hipFree(dcsrRowPtrB));
150    HIP_CHECK(hipFree(dcsrColIndB));
151    HIP_CHECK(hipFree(dcsrValB));
152    HIP_CHECK(hipFree(dcsrRowPtrC));
153    HIP_CHECK(hipFree(dcsrColIndC));
154    HIP_CHECK(hipFree(dcsrValC));
155
156    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrA));
157    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrB));
158    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrC));
159    HIPSPARSE_CHECK(hipsparseDestroy(handle));
160
161    return 0;
162}

hipsparseXcsrgemm2_bufferSizeExt()#

hipsparseStatus_t hipsparseScsrgemm2_bufferSizeExt(hipsparseHandle_t handle, int m, int n, int k, const float *alpha, const hipsparseMatDescr_t descrA, int nnzA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const int *csrRowPtrB, const int *csrColIndB, const float *beta, const hipsparseMatDescr_t descrD, int nnzD, const int *csrRowPtrD, const int *csrColIndD, csrgemm2Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDcsrgemm2_bufferSizeExt(hipsparseHandle_t handle, int m, int n, int k, const double *alpha, const hipsparseMatDescr_t descrA, int nnzA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const int *csrRowPtrB, const int *csrColIndB, const double *beta, const hipsparseMatDescr_t descrD, int nnzD, const int *csrRowPtrD, const int *csrColIndD, csrgemm2Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCcsrgemm2_bufferSizeExt(hipsparseHandle_t handle, int m, int n, int k, const hipComplex *alpha, const hipsparseMatDescr_t descrA, int nnzA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const int *csrRowPtrB, const int *csrColIndB, const hipComplex *beta, const hipsparseMatDescr_t descrD, int nnzD, const int *csrRowPtrD, const int *csrColIndD, csrgemm2Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZcsrgemm2_bufferSizeExt(hipsparseHandle_t handle, int m, int n, int k, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, int nnzA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const int *csrRowPtrB, const int *csrColIndB, const hipDoubleComplex *beta, const hipsparseMatDescr_t descrD, int nnzD, const int *csrRowPtrD, const int *csrColIndD, csrgemm2Info_t info, size_t *pBufferSizeInBytes)#

hipsparseXcsrgemm2_bufferSizeExt returns the size of the temporary storage buffer in bytes that is required by hipsparseXcsrgemm2Nnz() and hipsparseXcsrgemm2(). The temporary storage buffer must be allocated by the user.

Note

Please note, that for matrix products with more than 4096 non-zero entries per row, additional temporary storage buffer is allocated by the algorithm.

Note

Please note, that for matrix products with more than 8192 intermediate products per row, additional temporary storage buffer is allocated by the algorithm.

Note

Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

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

  • m[in] number of rows of the sparse CSR matrix \(op(A)\) and \(C\).

  • n[in] number of columns of the sparse CSR matrix \(op(B)\) and \(C\).

  • k[in] number of columns of the sparse CSR matrix \(op(A)\) and number of rows of the sparse CSR matrix \(op(B)\).

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

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

  • nnzA[in] number of non-zero entries of the sparse CSR matrix \(A\).

  • csrRowPtrA[in] array of m+1 elements ( \(op(A) == A\), k+1 otherwise) that point to the start of every row of the sparse CSR matrix \(op(A)\).

  • csrColIndA[in] array of nnzA elements containing the column indices of the sparse CSR matrix \(A\).

  • descrB[in] descriptor of the sparse CSR matrix \(B\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzB[in] number of non-zero entries of the sparse CSR matrix \(B\).

  • csrRowPtrB[in] array of k+1 elements ( \(op(B) == B\), m+1 otherwise) that point to the start of every row of the sparse CSR matrix \(op(B)\).

  • csrColIndB[in] array of nnzB elements containing the column indices of the sparse CSR matrix \(B\).

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

  • descrD[in] descriptor of the sparse CSR matrix \(D\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzD[in] number of non-zero entries of the sparse CSR matrix \(D\).

  • csrRowPtrD[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(D\).

  • csrColIndD[in] array of nnzD elements containing the column indices of the sparse CSR matrix \(D\).

  • info[inout] structure that holds meta data for the sparse CSR matrix \(C\).

  • pBufferSizeInBytes[out] number of bytes of the temporary storage buffer required by hipsparseXcsrgemm2Nnz(), hipsparseScsrgemm2(), hipsparseDcsrgemm2(), hipsparseCcsrgemm2() and hipsparseZcsrgemm2().

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, n, k, nnzA, nnzB, nnz_D, alpha, beta, descrA, csrRowPtrA, csrColIndA, descrB, csrRowPtrB, csrColIndB, descrD, csrRowPtrD, csrColIndD, info or pBufferSizeInBytes is invalid.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXcsrgemm2Nnz()#

hipsparseStatus_t hipsparseXcsrgemm2Nnz(hipsparseHandle_t handle, int m, int n, int k, const hipsparseMatDescr_t descrA, int nnzA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const int *csrRowPtrB, const int *csrColIndB, const hipsparseMatDescr_t descrD, int nnzD, const int *csrRowPtrD, const int *csrColIndD, const hipsparseMatDescr_t descrC, int *csrRowPtrC, int *nnzTotalDevHostPtr, const csrgemm2Info_t info, void *pBuffer)#

hipsparseXcsrgemm2Nnz computes the total CSR non-zero elements and the CSR row offsets, that point to the start of every row of the sparse CSR matrix, of the resulting multiplied matrix \(C\). It is assumed that csrRowPtrC has been allocated with size m+1. The required buffer size can be obtained by hipsparseXcsrgemm2_bufferSizeExt(). The desired index base in the output CSR matrix \(C\) is set in the hipsparseMatDescr_t descrC. See hipsparseSetMatIndexBase().

Note

As indicated, nnzTotalDevHostPtr can point either to host or device memory. This is controlled by setting the pointer mode. See hipsparseSetPointerMode().

Note

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

Note

Please note, that for matrix products with more than 8192 intermediate products per row, additional temporary storage buffer is allocated by the algorithm.

Note

Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

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

  • m[in] number of rows of the sparse CSR matrix \(op(A)\) and \(C\).

  • n[in] number of columns of the sparse CSR matrix \(op(B)\) and \(C\).

  • k[in] number of columns of the sparse CSR matrix \(op(A)\) and number of rows of the sparse CSR matrix \(op(B)\).

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

  • nnzA[in] number of non-zero entries of the sparse CSR matrix \(A\).

  • csrRowPtrA[in] array of m+1 elements ( \(op(A) == A\), k+1 otherwise) that point to the start of every row of the sparse CSR matrix \(op(A)\).

  • csrColIndA[in] array of nnzA elements containing the column indices of the sparse CSR matrix \(A\).

  • descrB[in] descriptor of the sparse CSR matrix \(B\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzB[in] number of non-zero entries of the sparse CSR matrix \(B\).

  • csrRowPtrB[in] array of k+1 elements ( \(op(B) == B\), m+1 otherwise) that point to the start of every row of the sparse CSR matrix \(op(B)\).

  • csrColIndB[in] array of nnzB elements containing the column indices of the sparse CSR matrix \(B\).

  • descrD[in] descriptor of the sparse CSR matrix \(D\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzD[in] number of non-zero entries of the sparse CSR matrix \(D\).

  • csrRowPtrD[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(D\).

  • csrColIndD[in] array of nnzD elements containing the column indices of the sparse CSR matrix \(D\).

  • descrC[in] descriptor of the sparse CSR matrix \(C\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • csrRowPtrC[out] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(C\).

  • nnzTotalDevHostPtr[out] pointer to the number of non-zero entries of the sparse CSR matrix \(C\).

  • info[in] structure that holds meta data for the sparse CSR matrix \(C\).

  • pBuffer[in] temporary storage buffer allocated by the user, size is returned by hipsparseScsrgemm2_bufferSizeExt(), hipsparseDcsrgemm2_bufferSizeExt(), hipsparseZcsrgemm2_bufferSizeExt() or hipsparseZcsrgemm2_bufferSizeExt().

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, n, k, nnzA, nnzB, nnzD, descrA, csrRowPtrA, csrColIndA, descrB, csrRowPtrB, csrColIndB, descrD, csrRowPtrD, csrColIndD, descrC, csrRowPtrC, nnzTotalDevHostPtr, info or pBuffer is invalid.

  • HIPSPARSE_STATUS_ALLOC_FAILED – additional buffer for long rows could not be allocated.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXcsrgemm2()#

hipsparseStatus_t hipsparseScsrgemm2(hipsparseHandle_t handle, int m, int n, int k, const float *alpha, const hipsparseMatDescr_t descrA, int nnzA, const float *csrValA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const float *csrValB, const int *csrRowPtrB, const int *csrColIndB, const float *beta, const hipsparseMatDescr_t descrD, int nnzD, const float *csrValD, const int *csrRowPtrD, const int *csrColIndD, const hipsparseMatDescr_t descrC, float *csrValC, const int *csrRowPtrC, int *csrColIndC, const csrgemm2Info_t info, void *pBuffer)#
hipsparseStatus_t hipsparseDcsrgemm2(hipsparseHandle_t handle, int m, int n, int k, const double *alpha, const hipsparseMatDescr_t descrA, int nnzA, const double *csrValA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const double *csrValB, const int *csrRowPtrB, const int *csrColIndB, const double *beta, const hipsparseMatDescr_t descrD, int nnzD, const double *csrValD, const int *csrRowPtrD, const int *csrColIndD, const hipsparseMatDescr_t descrC, double *csrValC, const int *csrRowPtrC, int *csrColIndC, const csrgemm2Info_t info, void *pBuffer)#
hipsparseStatus_t hipsparseCcsrgemm2(hipsparseHandle_t handle, int m, int n, int k, const hipComplex *alpha, const hipsparseMatDescr_t descrA, int nnzA, const hipComplex *csrValA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const hipComplex *csrValB, const int *csrRowPtrB, const int *csrColIndB, const hipComplex *beta, const hipsparseMatDescr_t descrD, int nnzD, const hipComplex *csrValD, const int *csrRowPtrD, const int *csrColIndD, const hipsparseMatDescr_t descrC, hipComplex *csrValC, const int *csrRowPtrC, int *csrColIndC, const csrgemm2Info_t info, void *pBuffer)#
hipsparseStatus_t hipsparseZcsrgemm2(hipsparseHandle_t handle, int m, int n, int k, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, int nnzA, const hipDoubleComplex *csrValA, const int *csrRowPtrA, const int *csrColIndA, const hipsparseMatDescr_t descrB, int nnzB, const hipDoubleComplex *csrValB, const int *csrRowPtrB, const int *csrColIndB, const hipDoubleComplex *beta, const hipsparseMatDescr_t descrD, int nnzD, const hipDoubleComplex *csrValD, const int *csrRowPtrD, const int *csrColIndD, const hipsparseMatDescr_t descrC, hipDoubleComplex *csrValC, const int *csrRowPtrC, int *csrColIndC, const csrgemm2Info_t info, void *pBuffer)#

Sparse matrix sparse matrix multiplication using CSR storage format.

hipsparseXcsrgemm2 multiplies the scalar \(\alpha\) with the sparse \(m \times k\) matrix \(A\), defined in CSR storage format, and the sparse \(k \times n\) matrix \(B\), defined in CSR storage format, and adds the result to the sparse \(m \times n\) matrix \(D\) that is multiplied by \(\beta\). The final result is stored in the sparse \(m \times n\) matrix \(C\), defined in CSR storage format, such that

\[ C := \alpha \cdot A \cdot B + \beta \cdot D \]

This computation involves a multi step process. First the user must call hipsparseXcsrgemm2_bufferSizeExt() in order to determine the required user allocated temporary buffer size. The user then allocates this buffer and also allocates csrRowPtrC to have size m+1. Both the temporary storage buffer and csrRowPtrC array are then passed to hipsparseXcsrgemm2Nnz which fills in the csrRowPtrC array as well as computes the total number of nonzeros in C, nnzC. The user then allocates both arrays csrColIndC and csrValC to have size nnzC and calls hipsparseXcsrgemm2 to complete the computation. The desired index base in the output CSR matrix C is set in the hipsparseMatDescr_t descrC. See hipsparseSetMatIndexBase().

Note

If \(\alpha == 0\), then \(C = \beta \cdot D\) will be computed.

Note

If \(\beta == 0\), then \(C = \alpha \cdot A \cdot B\) will be computed.

Note

\(\alpha == beta == 0\) is invalid.

Note

Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

Note

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

Note

Please note, that for matrix products with more than 4096 non-zero entries per row, additional temporary storage buffer is allocated by the algorithm.

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

  • m[in] number of rows of the sparse CSR matrix \(op(A)\) and \(C\).

  • n[in] number of columns of the sparse CSR matrix \(op(B)\) and \(C\).

  • k[in] number of columns of the sparse CSR matrix \(op(A)\) and number of rows of the sparse CSR matrix \(op(B)\).

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

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

  • nnzA[in] number of non-zero entries of the sparse CSR matrix \(A\).

  • csrValA[in] array of nnzA elements of the sparse CSR matrix \(A\).

  • csrRowPtrA[in] array of m+1 elements ( \(op(A) == A\), k+1 otherwise) that point to the start of every row of the sparse CSR matrix \(op(A)\).

  • csrColIndA[in] array of nnzA elements containing the column indices of the sparse CSR matrix \(A\).

  • descrB[in] descriptor of the sparse CSR matrix \(B\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzB[in] number of non-zero entries of the sparse CSR matrix \(B\).

  • csrValB[in] array of nnzB elements of the sparse CSR matrix \(B\).

  • csrRowPtrB[in] array of k+1 elements ( \(op(B) == B\), m+1 otherwise) that point to the start of every row of the sparse CSR matrix \(op(B)\).

  • csrColIndB[in] array of nnzB elements containing the column indices of the sparse CSR matrix \(B\).

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

  • descrD[in] descriptor of the sparse CSR matrix \(D\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • nnzD[in] number of non-zero entries of the sparse CSR matrix \(D\).

  • csrValD[in] array of nnzD elements of the sparse CSR matrix \(D\).

  • csrRowPtrD[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(D\).

  • csrColIndD[in] array of nnzD elements containing the column indices of the sparse CSR matrix \(D\).

  • descrC[in] descriptor of the sparse CSR matrix \(C\). Currently, only HIPSPARSE_MATRIX_TYPE_GENERAL is supported.

  • csrValC[out] array of nnzC elements of the sparse CSR matrix \(C\).

  • csrRowPtrC[in] array of m+1 elements that point to the start of every row of the sparse CSR matrix \(C\).

  • csrColIndC[out] array of nnzC elements containing the column indices of the sparse CSR matrix \(C\).

  • info[in] structure that holds meta data for the sparse CSR matrix \(C\).

  • pBuffer[in] temporary storage buffer allocated by the user, size is returned by hipsparseScsrgemm2_bufferSizeExt(), hipsparseDcsrgemm2_bufferSizeExt(), hipsparseCcsrgemm2_bufferSizeExt() or hipsparseZcsrgemm2_bufferSizeExt().

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, n, k, nnzA, nnzB, nnzD, alpha, beta, descrA, csrValA, csrRowPtrA, csrColIndA, descrB, csrValB, csrRowPtrB, csrColIndB, descrD, csrValD, csrRowPtrD, csrColIndD, csrValC, csrRowPtrC, csrColIndC, info or pBuffer is invalid.

  • HIPSPARSE_STATUS_ALLOC_FAILED – additional buffer for long rows could not be allocated.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

  1int main(int argc, char* argv[])
  2{
  3    int m    = 4;
  4    int k    = 3;
  5    int n    = 2;
  6    int nnzA = 7;
  7    int nnzB = 3;
  8    int nnzD = 6;
  9
 10    float alpha{1.0f};
 11    float beta{1.0f};
 12
 13    // A, B, and C are mxk, kxn, and m×n
 14
 15    // A
 16    // 1 0 0
 17    // 3 4 0
 18    // 5 6 7
 19    // 0 0 9
 20    std::vector<int>   hcsrRowPtrA = {0, 1, 3, 6, 7};
 21    std::vector<int>   hcsrColIndA = {0, 0, 1, 0, 1, 2, 2};
 22    std::vector<float> hcsrValA    = {1.0f, 3.0f, 4.0f, 5.0f, 6.0f, 7.0f, 9.0f};
 23
 24    // B
 25    // 0 1
 26    // 1 0
 27    // 0 1
 28    std::vector<int>   hcsrRowPtrB = {0, 1, 2, 3};
 29    std::vector<int>   hcsrColIndB = {1, 0, 1};
 30    std::vector<float> hcsrValB    = {1.0f, 1.0f, 1.0f};
 31
 32    // D
 33    // 0 1
 34    // 2 3
 35    // 4 5
 36    // 0 6
 37    std::vector<int>   hcsrRowPtrD = {0, 1, 3, 5, 6};
 38    std::vector<int>   hcsrColIndD = {1, 0, 1, 0, 1, 1};
 39    std::vector<float> hcsrValD    = {1.0f, 2.0f, 3.0f, 4.0f, 5.0f, 6.0f};
 40
 41    // Device memory management: Allocate and copy A, B
 42    int*   dcsrRowPtrA;
 43    int*   dcsrColIndA;
 44    float* dcsrValA;
 45    int*   dcsrRowPtrB;
 46    int*   dcsrColIndB;
 47    float* dcsrValB;
 48    int*   dcsrRowPtrD;
 49    int*   dcsrColIndD;
 50    float* dcsrValD;
 51    int*   dcsrRowPtrC;
 52    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrA, (m + 1) * sizeof(int)));
 53    HIP_CHECK(hipMalloc((void**)&dcsrColIndA, nnzA * sizeof(int)));
 54    HIP_CHECK(hipMalloc((void**)&dcsrValA, nnzA * sizeof(float)));
 55    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrB, (k + 1) * sizeof(int)));
 56    HIP_CHECK(hipMalloc((void**)&dcsrColIndB, nnzB * sizeof(int)));
 57    HIP_CHECK(hipMalloc((void**)&dcsrValB, nnzB * sizeof(float)));
 58    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrD, (m + 1) * sizeof(int)));
 59    HIP_CHECK(hipMalloc((void**)&dcsrColIndD, nnzD * sizeof(int)));
 60    HIP_CHECK(hipMalloc((void**)&dcsrValD, nnzD * sizeof(float)));
 61    HIP_CHECK(hipMalloc((void**)&dcsrRowPtrC, (m + 1) * sizeof(int)));
 62
 63    HIP_CHECK(
 64        hipMemcpy(dcsrRowPtrA, hcsrRowPtrA.data(), (m + 1) * sizeof(int), hipMemcpyHostToDevice));
 65    HIP_CHECK(
 66        hipMemcpy(dcsrColIndA, hcsrColIndA.data(), nnzA * sizeof(int), hipMemcpyHostToDevice));
 67    HIP_CHECK(hipMemcpy(dcsrValA, hcsrValA.data(), nnzA * sizeof(float), hipMemcpyHostToDevice));
 68    HIP_CHECK(
 69        hipMemcpy(dcsrRowPtrB, hcsrRowPtrB.data(), (k + 1) * sizeof(int), hipMemcpyHostToDevice));
 70    HIP_CHECK(
 71        hipMemcpy(dcsrColIndB, hcsrColIndB.data(), nnzB * sizeof(int), hipMemcpyHostToDevice));
 72    HIP_CHECK(hipMemcpy(dcsrValB, hcsrValB.data(), nnzB * sizeof(float), hipMemcpyHostToDevice));
 73    HIP_CHECK(
 74        hipMemcpy(dcsrRowPtrD, hcsrRowPtrD.data(), (m + 1) * sizeof(int), hipMemcpyHostToDevice));
 75    HIP_CHECK(
 76        hipMemcpy(dcsrColIndD, hcsrColIndD.data(), nnzD * sizeof(int), hipMemcpyHostToDevice));
 77    HIP_CHECK(hipMemcpy(dcsrValD, hcsrValD.data(), nnzD * sizeof(float), hipMemcpyHostToDevice));
 78
 79    hipsparseHandle_t handle;
 80    HIPSPARSE_CHECK(hipsparseCreate(&handle));
 81
 82    hipsparseMatDescr_t descrA;
 83    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrA));
 84
 85    hipsparseMatDescr_t descrB;
 86    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrB));
 87
 88    hipsparseMatDescr_t descrC;
 89    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrC));
 90
 91    hipsparseMatDescr_t descrD;
 92    HIPSPARSE_CHECK(hipsparseCreateMatDescr(&descrD));
 93
 94    csrgemm2Info_t info;
 95    HIPSPARSE_CHECK(hipsparseCreateCsrgemm2Info(&info));
 96
 97    size_t bufferSize;
 98    HIPSPARSE_CHECK(hipsparseScsrgemm2_bufferSizeExt(handle,
 99                                                     m,
100                                                     n,
101                                                     k,
102                                                     &alpha,
103                                                     descrA,
104                                                     nnzA,
105                                                     dcsrRowPtrA,
106                                                     dcsrColIndA,
107                                                     descrB,
108                                                     nnzB,
109                                                     dcsrRowPtrB,
110                                                     dcsrColIndB,
111                                                     &beta,
112                                                     descrD,
113                                                     nnzD,
114                                                     dcsrRowPtrD,
115                                                     dcsrColIndD,
116                                                     info,
117                                                     &bufferSize));
118
119    void* dbuffer = nullptr;
120    HIP_CHECK(hipMalloc((void**)&dbuffer, bufferSize));
121
122    int nnzC;
123    HIPSPARSE_CHECK(hipsparseXcsrgemm2Nnz(handle,
124                                          m,
125                                          n,
126                                          k,
127                                          descrA,
128                                          nnzA,
129                                          dcsrRowPtrA,
130                                          dcsrColIndA,
131                                          descrB,
132                                          nnzB,
133                                          dcsrRowPtrB,
134                                          dcsrColIndB,
135                                          descrD,
136                                          nnzD,
137                                          dcsrRowPtrD,
138                                          dcsrColIndD,
139                                          descrC,
140                                          dcsrRowPtrC,
141                                          &nnzC,
142                                          info,
143                                          dbuffer));
144
145    int*   dcsrColIndC = nullptr;
146    float* dcsrValC    = nullptr;
147    HIP_CHECK(hipMalloc((void**)&dcsrColIndC, sizeof(int) * nnzC));
148    HIP_CHECK(hipMalloc((void**)&dcsrValC, sizeof(float) * nnzC));
149
150    HIPSPARSE_CHECK(hipsparseScsrgemm2(handle,
151                                       m,
152                                       n,
153                                       k,
154                                       &alpha,
155                                       descrA,
156                                       nnzA,
157                                       dcsrValA,
158                                       dcsrRowPtrA,
159                                       dcsrColIndA,
160                                       descrB,
161                                       nnzB,
162                                       dcsrValB,
163                                       dcsrRowPtrB,
164                                       dcsrColIndB,
165                                       &beta,
166                                       descrD,
167                                       nnzD,
168                                       dcsrValD,
169                                       dcsrRowPtrD,
170                                       dcsrColIndD,
171                                       descrC,
172                                       dcsrValC,
173                                       dcsrRowPtrC,
174                                       dcsrColIndC,
175                                       info,
176                                       dbuffer));
177
178    std::vector<int>   hcsrRowPtrC(m + 1);
179    std::vector<int>   hcsrColIndC(nnzC);
180    std::vector<float> hcsrValC(nnzC);
181
182    // Copy back to the host
183    HIP_CHECK(
184        hipMemcpy(hcsrRowPtrC.data(), dcsrRowPtrC, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
185    HIP_CHECK(
186        hipMemcpy(hcsrColIndC.data(), dcsrColIndC, sizeof(int) * nnzC, hipMemcpyDeviceToHost));
187    HIP_CHECK(hipMemcpy(hcsrValC.data(), dcsrValC, sizeof(float) * nnzC, hipMemcpyDeviceToHost));
188
189    std::cout << "C" << std::endl;
190    for(int i = 0; i < m; i++)
191    {
192        int start = hcsrRowPtrC[i];
193        int end   = hcsrRowPtrC[i + 1];
194
195        std::vector<float> temp(n, 0.0f);
196        for(int j = start; j < end; j++)
197        {
198            temp[hcsrColIndC[j]] = hcsrValC[j];
199        }
200
201        for(int j = 0; j < n; j++)
202        {
203            std::cout << temp[j] << " ";
204        }
205        std::cout << std::endl;
206    }
207    std::cout << std::endl;
208
209    HIP_CHECK(hipFree(dcsrRowPtrA));
210    HIP_CHECK(hipFree(dcsrColIndA));
211    HIP_CHECK(hipFree(dcsrValA));
212    HIP_CHECK(hipFree(dcsrRowPtrB));
213    HIP_CHECK(hipFree(dcsrColIndB));
214    HIP_CHECK(hipFree(dcsrValB));
215    HIP_CHECK(hipFree(dcsrRowPtrC));
216    HIP_CHECK(hipFree(dcsrColIndC));
217    HIP_CHECK(hipFree(dcsrValC));
218    HIP_CHECK(hipFree(dcsrRowPtrD));
219    HIP_CHECK(hipFree(dcsrColIndD));
220    HIP_CHECK(hipFree(dcsrValD));
221
222    HIP_CHECK(hipFree(dbuffer));
223
224    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrA));
225    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrB));
226    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrC));
227    HIPSPARSE_CHECK(hipsparseDestroyMatDescr(descrD));
228    HIPSPARSE_CHECK(hipsparseDestroyCsrgemm2Info(info));
229
230    HIPSPARSE_CHECK(hipsparseDestroy(handle));
231
232    return 0;
233}