Sparse generic functions#
This module contains all sparse generic routines.
The sparse generic routines describe operations that manipulate sparse matrices.
hipsparseAxpby()#
-
hipsparseStatus_t hipsparseAxpby(hipsparseHandle_t handle, const void *alpha, hipsparseConstSpVecDescr_t vecX, const void *beta, hipsparseDnVecDescr_t vecY)#
Scale a sparse vector and add it to a scaled dense vector.
hipsparseAxpbymultiplies the sparse vector \(x\) with scalar \(\alpha\) and adds the result to the dense vector \(y\) that is multiplied with scalar \(\beta\), such that\[ y := \alpha \cdot x + \beta \cdot y \]for(i = 0; i < size; ++i) { y[i] = beta * y[i] } for(i = 0; i < nnz; ++i) { y[xInd[i]] += alpha * xVal[i] }
hipsparseAxpbysupports the following precision data types for the sparse and dense vectors \(x\) and \(y\) and compute types for the scalars \(\alpha\) and \(\beta\).- Uniform Precisions:
X / Y / compute_type
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
- Mixed precisions:
X / Y
compute_type
HIP_R_16F
HIP_R_32F
HIP_R_16BF
HIP_R_32F
- Parameters:
handle – [in] handle to the hipsparse library context queue.
alpha – [in] scalar \(\alpha\).
vecX – [in] sparse vector descriptor.
beta – [in] scalar \(\beta\).
vecY – [inout] dense vector descriptor.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_NOT_INITIALIZED –
handleis not initialized.HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,vecX,betaorvecYis nullptr, or the vector sizes or data types are incompatible.
1int main(int argc, char* argv[])
2{
3 // Number of non-zeros of the sparse vector
4 int nnz = 3;
5
6 // Size of sparse and dense vector
7 int size = 9;
8
9 // Sparse index vector
10 std::vector<int> hxInd = {0, 3, 5};
11
12 // Sparse value vector
13 std::vector<float> hxVal = {1.0f, 2.0f, 3.0f};
14
15 // Dense vector
16 std::vector<float> hy = {1.0f, 2.0f, 3.0f, 4.0f, 5.0f, 6.0f, 7.0f, 8.0f, 9.0f};
17
18 // Scalar alpha
19 float alpha = 3.7f;
20
21 // Scalar beta
22 float beta = 1.2f;
23
24 // Offload data to device
25 int* dxInd;
26 float* dxVal;
27 float* dy;
28 HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
29 HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
30 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
31
32 HIP_CHECK(hipMemcpy(dxInd, hxInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
33 HIP_CHECK(hipMemcpy(dxVal, hxVal.data(), sizeof(float) * nnz, hipMemcpyHostToDevice));
34 HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(float) * size, hipMemcpyHostToDevice));
35
36 hipsparseHandle_t handle;
37 HIPSPARSE_CHECK(hipsparseCreate(&handle));
38
39 // Create sparse vector X
40 hipsparseSpVecDescr_t vecX;
41 HIPSPARSE_CHECK(hipsparseCreateSpVec(
42 &vecX, size, nnz, dxInd, dxVal, HIPSPARSE_INDEX_32I, HIPSPARSE_INDEX_BASE_ZERO, HIP_R_32F));
43
44 // Create dense vector Y
45 hipsparseDnVecDescr_t vecY;
46 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, size, dy, HIP_R_32F));
47
48 // Call axpby to perform y = beta * y + alpha * x
49 HIPSPARSE_CHECK(hipsparseAxpby(handle, &alpha, vecX, &beta, vecY));
50
51 HIPSPARSE_CHECK(hipsparseDnVecGetValues(vecY, (void**)&dy));
52
53 // Copy result back to host
54 HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(float) * size, hipMemcpyDeviceToHost));
55
56 std::cout << "hy" << std::endl;
57 for(size_t i = 0; i < hy.size(); i++)
58 {
59 std::cout << hy[i] << " ";
60 }
61 std::cout << std::endl;
62
63 // Clear hipSPARSE
64 HIPSPARSE_CHECK(hipsparseDestroySpVec(vecX));
65 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
66 HIPSPARSE_CHECK(hipsparseDestroy(handle));
67
68 // Clear device memory
69 HIP_CHECK(hipFree(dxInd));
70 HIP_CHECK(hipFree(dxVal));
71 HIP_CHECK(hipFree(dy));
72
73 return 0;
74}
1int main(int argc, char* argv[])
2{
3 // Number of non-zeros of the sparse vector
4 int nnz = 3;
5
6 // Size of sparse and dense vector
7 int size = 9;
8
9 // Sparse index vector
10 int hxInd[] = {0, 3, 5};
11
12 // Sparse value vector
13 float hxVal[] = {1.0, 2.0, 3.0};
14
15 // Dense vector
16 float hy[] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0};
17
18 // Scalar alpha
19 float alpha = 3.7;
20
21 // Scalar beta
22 float beta = 1.2;
23
24 // Offload data to device
25 int* dxInd;
26 float* dxVal;
27 float* dy;
28 HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
29 HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
30 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
31
32 HIP_CHECK(hipMemcpy(dxInd, hxInd, sizeof(int) * nnz, hipMemcpyHostToDevice));
33 HIP_CHECK(hipMemcpy(dxVal, hxVal, sizeof(float) * nnz, hipMemcpyHostToDevice));
34 HIP_CHECK(hipMemcpy(dy, hy, sizeof(float) * size, hipMemcpyHostToDevice));
35
36 hipsparseHandle_t handle;
37 HIPSPARSE_CHECK(hipsparseCreate(&handle));
38
39 // Create sparse vector X
40 hipsparseSpVecDescr_t vecX;
41 HIPSPARSE_CHECK(hipsparseCreateSpVec(
42 &vecX, size, nnz, dxInd, dxVal, HIPSPARSE_INDEX_32I, HIPSPARSE_INDEX_BASE_ZERO, HIP_R_32F));
43
44 // Create dense vector Y
45 hipsparseDnVecDescr_t vecY;
46 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, size, dy, HIP_R_32F));
47
48 // Call axpby to perform y = beta * y + alpha * x
49 HIPSPARSE_CHECK(hipsparseAxpby(handle, &alpha, vecX, &beta, vecY));
50
51 HIPSPARSE_CHECK(hipsparseDnVecGetValues(vecY, (void**)&dy));
52
53 // Copy result back to host
54 HIP_CHECK(hipMemcpy(hy, dy, sizeof(float) * size, hipMemcpyDeviceToHost));
55
56 printf("hy\n");
57 for(size_t i = 0; i < size; i++)
58 {
59 printf("%f ", hy[i]);
60 }
61 printf("\n");
62
63 // Clear hipSPARSE
64 HIPSPARSE_CHECK(hipsparseDestroySpVec(vecX));
65 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
66 HIPSPARSE_CHECK(hipsparseDestroy(handle));
67
68 // Clear device memory
69 HIP_CHECK(hipFree(dxInd));
70 HIP_CHECK(hipFree(dxVal));
71 HIP_CHECK(hipFree(dy));
72
73 return 0;
74}
hipsparseGather()#
-
hipsparseStatus_t hipsparseGather(hipsparseHandle_t handle, hipsparseConstDnVecDescr_t vecY, hipsparseSpVecDescr_t vecX)#
Gather elements from a dense vector and store them into a sparse vector.
hipsparseGathergathers the elements from the dense vector \(y\) and stores them in the sparse vector \(x\).for(i = 0; i < nnz; ++i) { x_val[i] = y[x_ind[i]]; }
hipsparseGathersupports the following uniform precision data types for the sparse and dense vectors \(x\) and \(y\).- Uniform Precisions:
X / Y
HIP_R_8I
HIP_R_16F
HIP_R_16BF
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
- Parameters:
handle – [in] handle to the hipsparse library context queue.
vecY – [in] dense vector descriptor \(y\).
vecX – [out] sparse vector descriptor \(x\).
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_NOT_INITIALIZED –
handleis not initialized.HIPSPARSE_STATUS_INVALID_VALUE –
handle,vecXorvecYis nullptr, or the vector sizes or data types are incompatible.
1int main(int argc, char* argv[])
2{
3 // Number of non-zeros of the sparse vector
4 int nnz = 3;
5
6 // Size of sparse and dense vector
7 int size = 9;
8
9 // Sparse index vector
10 std::vector<int> hxInd = {0, 3, 5};
11
12 // Dense vector
13 std::vector<float> hy = {1.0f, 2.0f, 3.0f, 4.0f, 5.0f, 6.0f, 7.0f, 8.0f, 9.0f};
14
15 // Offload data to device
16 int* dxInd;
17 float* dxVal;
18 float* dy;
19 HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
20 HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
21 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
22
23 HIP_CHECK(hipMemcpy(dxInd, hxInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
24 HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(float) * size, hipMemcpyHostToDevice));
25
26 hipsparseHandle_t handle;
27 HIPSPARSE_CHECK(hipsparseCreate(&handle));
28
29 // Create sparse vector X
30 hipsparseSpVecDescr_t vecX;
31 HIPSPARSE_CHECK(hipsparseCreateSpVec(
32 &vecX, size, nnz, dxInd, dxVal, HIPSPARSE_INDEX_32I, HIPSPARSE_INDEX_BASE_ZERO, HIP_R_32F));
33
34 // Create dense vector Y
35 hipsparseDnVecDescr_t vecY;
36 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, size, dy, HIP_R_32F));
37
38 // Perform gather
39 HIPSPARSE_CHECK(hipsparseGather(handle, vecY, vecX));
40
41 HIPSPARSE_CHECK(hipsparseSpVecGetValues(vecX, (void**)&dxVal));
42
43 // Copy result back to host
44 std::vector<float> hxVal(nnz, 0.0f);
45 HIP_CHECK(hipMemcpy(hxVal.data(), dxVal, sizeof(float) * nnz, hipMemcpyDeviceToHost));
46
47 // Clear hipSPARSE
48 HIPSPARSE_CHECK(hipsparseDestroySpVec(vecX));
49 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
50 HIPSPARSE_CHECK(hipsparseDestroy(handle));
51
52 // Clear device memory
53 HIP_CHECK(hipFree(dxInd));
54 HIP_CHECK(hipFree(dxVal));
55 HIP_CHECK(hipFree(dy));
56
57 return 0;
58}
1int main(int argc, char* argv[])
2{
3 // Number of non-zeros of the sparse vector
4 int nnz = 3;
5
6 // Size of sparse and dense vector
7 int size = 9;
8
9 // Sparse index vector
10 int hxInd[] = {0, 3, 5};
11
12 // Dense vector
13 float hy[] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0};
14
15 // Offload data to device
16 int* dxInd;
17 float* dxVal;
18 float* dy;
19 HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
20 HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
21 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
22
23 HIP_CHECK(hipMemcpy(dxInd, hxInd, sizeof(int) * nnz, hipMemcpyHostToDevice));
24 HIP_CHECK(hipMemcpy(dy, hy, sizeof(float) * size, hipMemcpyHostToDevice));
25
26 hipsparseHandle_t handle;
27 HIPSPARSE_CHECK(hipsparseCreate(&handle));
28
29 // Create sparse vector X
30 hipsparseSpVecDescr_t vecX;
31 HIPSPARSE_CHECK(hipsparseCreateSpVec(
32 &vecX, size, nnz, dxInd, dxVal, HIPSPARSE_INDEX_32I, HIPSPARSE_INDEX_BASE_ZERO, HIP_R_32F));
33
34 // Create dense vector Y
35 hipsparseDnVecDescr_t vecY;
36 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, size, dy, HIP_R_32F));
37
38 // Perform gather
39 HIPSPARSE_CHECK(hipsparseGather(handle, vecY, vecX));
40
41 HIPSPARSE_CHECK(hipsparseSpVecGetValues(vecX, (void**)&dxVal));
42
43 // Copy result back to host
44 float* hxVal = (float*)malloc(nnz * sizeof(float));
45 HIP_CHECK(hipMemcpy(hxVal, dxVal, sizeof(float) * nnz, hipMemcpyDeviceToHost));
46
47 // Clear hipSPARSE
48 HIPSPARSE_CHECK(hipsparseDestroySpVec(vecX));
49 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
50 HIPSPARSE_CHECK(hipsparseDestroy(handle));
51
52 // Clear device memory
53 HIP_CHECK(hipFree(dxInd));
54 HIP_CHECK(hipFree(dxVal));
55 HIP_CHECK(hipFree(dy));
56
57 return 0;
58}
hipsparseScatter()#
-
hipsparseStatus_t hipsparseScatter(hipsparseHandle_t handle, hipsparseConstSpVecDescr_t vecX, hipsparseDnVecDescr_t vecY)#
Scatter elements from a sparse vector into a dense vector.
hipsparseScatterscatters the elements from the sparse vector \(x\) in the dense vector \(y\).for(i = 0; i < nnz; ++i) { y[x_ind[i]] = x_val[i]; }
hipsparseScattersupports the following uniform precision data types for the sparse and dense vectors \(x\) and \(y\).- Uniform Precisions:
X / Y
HIP_R_8I
HIP_R_16F
HIP_R_16BF
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
- Parameters:
handle – [in] handle to the hipsparse library context queue.
vecX – [in] sparse vector descriptor \(x\).
vecY – [out] dense vector descriptor \(y\).
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_NOT_INITIALIZED –
handleis not initialized.HIPSPARSE_STATUS_INVALID_VALUE –
handle,vecXorvecYis nullptr, or the vector sizes or data types are incompatible.
1int main(int argc, char* argv[])
2{
3 // Number of non-zeros of the sparse vector
4 int nnz = 3;
5
6 // Size of sparse and dense vector
7 int size = 9;
8
9 // Sparse index vector
10 std::vector<int> hxInd = {0, 3, 5};
11
12 // Sparse value vector
13 std::vector<float> hxVal = {1.0f, 2.0f, 3.0f};
14
15 // Dense vector
16 std::vector<float> hy = {1.0f, 2.0f, 3.0f, 4.0f, 5.0f, 6.0f, 7.0f, 8.0f, 9.0f};
17
18 // Offload data to device
19 int* dxInd;
20 float* dxVal;
21 float* dy;
22 HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
23 HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
24 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
25
26 HIP_CHECK(hipMemcpy(dxInd, hxInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
27 HIP_CHECK(hipMemcpy(dxVal, hxVal.data(), sizeof(float) * nnz, hipMemcpyHostToDevice));
28 HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(float) * size, hipMemcpyHostToDevice));
29
30 hipsparseHandle_t handle;
31 HIPSPARSE_CHECK(hipsparseCreate(&handle));
32
33 // Create sparse vector X
34 hipsparseSpVecDescr_t vecX;
35 HIPSPARSE_CHECK(hipsparseCreateSpVec(
36 &vecX, size, nnz, dxInd, dxVal, HIPSPARSE_INDEX_32I, HIPSPARSE_INDEX_BASE_ZERO, HIP_R_32F));
37
38 // Create dense vector Y
39 hipsparseDnVecDescr_t vecY;
40 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, size, dy, HIP_R_32F));
41
42 // Perform scatter
43 HIPSPARSE_CHECK(hipsparseScatter(handle, vecX, vecY));
44
45 HIPSPARSE_CHECK(hipsparseDnVecGetValues(vecY, (void**)&dy));
46
47 // Copy result back to host
48 HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(float) * size, hipMemcpyDeviceToHost));
49
50 // Clear hipSPARSE
51 HIPSPARSE_CHECK(hipsparseDestroySpVec(vecX));
52 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
53 HIPSPARSE_CHECK(hipsparseDestroy(handle));
54
55 // Clear device memory
56 HIP_CHECK(hipFree(dxInd));
57 HIP_CHECK(hipFree(dxVal));
58 HIP_CHECK(hipFree(dy));
59
60 return 0;
61}
1int main(int argc, char* argv[])
2{
3 // Number of non-zeros of the sparse vector
4 int nnz = 3;
5
6 // Size of sparse and dense vector
7 int size = 9;
8
9 // Sparse index vector
10 int hxInd[] = {0, 3, 5};
11
12 // Sparse value vector
13 float hxVal[] = {1.0, 2.0, 3.0};
14
15 // Dense vector
16 float hy[] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0};
17
18 // Offload data to device
19 int* dxInd;
20 float* dxVal;
21 float* dy;
22 HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
23 HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
24 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
25
26 HIP_CHECK(hipMemcpy(dxInd, hxInd, sizeof(int) * nnz, hipMemcpyHostToDevice));
27 HIP_CHECK(hipMemcpy(dxVal, hxVal, sizeof(float) * nnz, hipMemcpyHostToDevice));
28 HIP_CHECK(hipMemcpy(dy, hy, sizeof(float) * size, hipMemcpyHostToDevice));
29
30 hipsparseHandle_t handle;
31 HIPSPARSE_CHECK(hipsparseCreate(&handle));
32
33 // Create sparse vector X
34 hipsparseSpVecDescr_t vecX;
35 HIPSPARSE_CHECK(hipsparseCreateSpVec(
36 &vecX, size, nnz, dxInd, dxVal, HIPSPARSE_INDEX_32I, HIPSPARSE_INDEX_BASE_ZERO, HIP_R_32F));
37
38 // Create dense vector Y
39 hipsparseDnVecDescr_t vecY;
40 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, size, dy, HIP_R_32F));
41
42 // Perform scatter
43 HIPSPARSE_CHECK(hipsparseScatter(handle, vecX, vecY));
44
45 HIPSPARSE_CHECK(hipsparseDnVecGetValues(vecY, (void**)&dy));
46
47 // Copy result back to host
48 HIP_CHECK(hipMemcpy(hy, dy, sizeof(float) * size, hipMemcpyDeviceToHost));
49
50 // Clear hipSPARSE
51 HIPSPARSE_CHECK(hipsparseDestroySpVec(vecX));
52 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
53 HIPSPARSE_CHECK(hipsparseDestroy(handle));
54
55 // Clear device memory
56 HIP_CHECK(hipFree(dxInd));
57 HIP_CHECK(hipFree(dxVal));
58 HIP_CHECK(hipFree(dy));
59
60 return 0;
61}
hipsparseRot()#
-
hipsparseStatus_t hipsparseRot(hipsparseHandle_t handle, const void *c_coeff, const void *s_coeff, hipsparseSpVecDescr_t vecX, hipsparseDnVecDescr_t vecY)#
Apply Givens rotation to a dense and a sparse vector.
hipsparseRotapplies the Givens rotation matrix \(G\) to the sparse vector \(x\) and the dense vector \(y\), where\[\begin{split} G = \begin{pmatrix} c & s \\ -s & c \end{pmatrix} \end{split}\]for(i = 0; i < nnz; ++i) { x_tmp = x_val[i]; y_tmp = y[x_ind[i]]; x_val[i] = c * x_tmp + s * y_tmp; y[x_ind[i]] = c * y_tmp - s * x_tmp; }
hipsparseRotsupports the following uniform precision data types for the sparse and dense vectors \(x\) and \(y\) and compute types for the scalars \(c\) and \(s\).- 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.
- Uniform Precisions:
X / Y / compute_type
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
- Parameters:
handle – [in] handle to the hipsparse library context queue.
c_coeff – [in] pointer to the cosine element of \(G\), can be on host or device.
s_coeff – [in] pointer to the sine element of \(G\), can be on host or device.
vecX – [inout] sparse vector descriptor \(x\).
vecY – [inout] dense vector descriptor \(y\).
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_NOT_INITIALIZED –
handleis not initialized.HIPSPARSE_STATUS_INVALID_VALUE –
handle,c_coeff,s_coeff,vecXorvecYis nullptr, or the vector sizes or data types are incompatible.
1int main(int argc, char* argv[])
2{
3 // Number of non-zeros of the sparse vector
4 int nnz = 3;
5
6 // Size of sparse and dense vector
7 int size = 9;
8
9 // Sparse index vector
10 std::vector<int> hxInd = {0, 3, 5};
11
12 // Sparse value vector
13 std::vector<float> hxVal = {1.0f, 2.0f, 3.0f};
14
15 // Dense vector
16 std::vector<float> hy = {1.0f, 2.0f, 3.0f, 4.0f, 5.0f, 6.0f, 7.0f, 8.0f, 9.0f};
17
18 // Scalar c
19 float c = 3.7f;
20
21 // Scalar s
22 float s = 1.2f;
23
24 // Offload data to device
25 int* dxInd;
26 float* dxVal;
27 float* dy;
28 HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
29 HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
30 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
31
32 HIP_CHECK(hipMemcpy(dxInd, hxInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
33 HIP_CHECK(hipMemcpy(dxVal, hxVal.data(), sizeof(float) * nnz, hipMemcpyHostToDevice));
34 HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(float) * size, hipMemcpyHostToDevice));
35
36 hipsparseHandle_t handle;
37 HIPSPARSE_CHECK(hipsparseCreate(&handle));
38
39 // Create sparse vector X
40 hipsparseSpVecDescr_t vecX;
41 HIPSPARSE_CHECK(hipsparseCreateSpVec(
42 &vecX, size, nnz, dxInd, dxVal, HIPSPARSE_INDEX_32I, HIPSPARSE_INDEX_BASE_ZERO, HIP_R_32F));
43
44 // Create dense vector Y
45 hipsparseDnVecDescr_t vecY;
46 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, size, dy, HIP_R_32F));
47
48 // Call rot
49 HIPSPARSE_CHECK(hipsparseRot(handle, (void*)&c, (void*)&s, vecX, vecY));
50
51 HIPSPARSE_CHECK(hipsparseSpVecGetValues(vecX, (void**)&dxVal));
52 HIPSPARSE_CHECK(hipsparseDnVecGetValues(vecY, (void**)&dy));
53
54 // Copy result back to host
55 HIP_CHECK(hipMemcpy(hxVal.data(), dxVal, sizeof(float) * nnz, hipMemcpyDeviceToHost));
56 HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(float) * size, hipMemcpyDeviceToHost));
57
58 // Clear hipSPARSE
59 HIPSPARSE_CHECK(hipsparseDestroySpVec(vecX));
60 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
61 HIPSPARSE_CHECK(hipsparseDestroy(handle));
62
63 // Clear device memory
64 HIP_CHECK(hipFree(dxInd));
65 HIP_CHECK(hipFree(dxVal));
66 HIP_CHECK(hipFree(dy));
67
68 return 0;
69}
1int main(int argc, char* argv[])
2{
3 // Number of non-zeros of the sparse vector
4 int nnz = 3;
5
6 // Size of sparse and dense vector
7 int size = 9;
8
9 // Sparse index vector
10 int hxInd[] = {0, 3, 5};
11
12 // Sparse value vector
13 float hxVal[] = {1.0, 2.0, 3.0};
14
15 // Dense vector
16 float hy[] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0};
17
18 // Scalar c
19 float c = 3.7;
20
21 // Scalar s
22 float s = 1.2;
23
24 // Offload data to device
25 int* dxInd;
26 float* dxVal;
27 float* dy;
28 HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
29 HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
30 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
31
32 HIP_CHECK(hipMemcpy(dxInd, hxInd, sizeof(int) * nnz, hipMemcpyHostToDevice));
33 HIP_CHECK(hipMemcpy(dxVal, hxVal, sizeof(float) * nnz, hipMemcpyHostToDevice));
34 HIP_CHECK(hipMemcpy(dy, hy, sizeof(float) * size, hipMemcpyHostToDevice));
35
36 hipsparseHandle_t handle;
37 HIPSPARSE_CHECK(hipsparseCreate(&handle));
38
39 // Create sparse vector X
40 hipsparseSpVecDescr_t vecX;
41 HIPSPARSE_CHECK(hipsparseCreateSpVec(
42 &vecX, size, nnz, dxInd, dxVal, HIPSPARSE_INDEX_32I, HIPSPARSE_INDEX_BASE_ZERO, HIP_R_32F));
43
44 // Create dense vector Y
45 hipsparseDnVecDescr_t vecY;
46 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, size, dy, HIP_R_32F));
47
48 // Call rot
49 HIPSPARSE_CHECK(hipsparseRot(handle, (void*)&c, (void*)&s, vecX, vecY));
50
51 HIPSPARSE_CHECK(hipsparseSpVecGetValues(vecX, (void**)&dxVal));
52 HIPSPARSE_CHECK(hipsparseDnVecGetValues(vecY, (void**)&dy));
53
54 // Copy result back to host
55 HIP_CHECK(hipMemcpy(hxVal, dxVal, sizeof(float) * nnz, hipMemcpyDeviceToHost));
56 HIP_CHECK(hipMemcpy(hy, dy, sizeof(float) * size, hipMemcpyDeviceToHost));
57
58 // Clear hipSPARSE
59 HIPSPARSE_CHECK(hipsparseDestroySpVec(vecX));
60 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
61 HIPSPARSE_CHECK(hipsparseDestroy(handle));
62
63 // Clear device memory
64 HIP_CHECK(hipFree(dxInd));
65 HIP_CHECK(hipFree(dxVal));
66 HIP_CHECK(hipFree(dy));
67
68 return 0;
69}
hipsparseSparseToDense_bufferSize()#
-
hipsparseStatus_t hipsparseSparseToDense_bufferSize(hipsparseHandle_t handle, hipsparseConstSpMatDescr_t matA, hipsparseDnMatDescr_t matB, hipsparseSparseToDenseAlg_t alg, size_t *pBufferSizeInBytes)#
hipsparseSparseToDense_bufferSizecomputes the required user allocated buffer size needed when converting a sparse matrix to a dense matrix. This routine currently accepts the sparse matrix descriptormatAin CSR, CSC, or COO format. This routine is used to determine the size of the buffer needed in hipsparseSparseToDense.hipsparseSparseToDense_bufferSizesupports different data types for the sparse and dense matrices. See hipsparseSparseToDense for a complete listing of all the data types available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
matA – [in] sparse matrix descriptor.
matB – [in] dense matrix descriptor.
alg – [in] algorithm for the sparse to dense computation.
pBufferSizeInBytes – [out] number of bytes of the temporary storage buffer.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,matA,matB, orpBufferSizeInBytespointer is invalid.
hipsparseSparseToDense()#
-
hipsparseStatus_t hipsparseSparseToDense(hipsparseHandle_t handle, hipsparseConstSpMatDescr_t matA, hipsparseDnMatDescr_t matB, hipsparseSparseToDenseAlg_t alg, void *externalBuffer)#
Sparse matrix to dense matrix conversion.
hipsparseSparseToDenseconverts a sparse matrix to a dense matrix. This routine currently accepts the sparse matrix descriptormatAin CSR, CSC, or COO format. This routine takes a user allocated buffer whose size must first be computed by calling hipsparseSparseToDense_bufferSizeThe conversion of a sparse matrix into a dense one involves two steps. First, the user creates the sparse and dense matrix descriptors and calls hipsparseSparseToDense_bufferSize to determine the size of the required temporary storage buffer. The user then allocates this buffer and passes it to hipsparseSparseToDense in order to complete the conversion. Once the conversion is complete, the user is free to deallocate the storage buffer. See full example below for details.
hipsparseSparseToDensesupports the following uniform precision data types for the sparse and dense matrices \(A\) and \(B\):- Uniform Precisions:
A / B
HIP_R_16F
HIP_R_16BF
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
Note
Currently only the sparse matrix formats CSR, CSC, and COO are supported when converting a sparse matrix to a dense matrix.
- Parameters:
handle – [in] handle to the hipsparse library context queue.
matA – [in] sparse matrix descriptor.
matB – [in] dense matrix descriptor.
alg – [in] algorithm for the sparse to dense computation.
externalBuffer – [in] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,matA,matB, orexternalBufferpointer is invalid.
1int main(int argc, char* argv[])
2{
3 // 1 0 0 0
4 // A = 4 2 0 4
5 // 0 3 7 0
6 // 9 0 0 1
7 int m = 4;
8 int n = 4;
9 int nnz = 8;
10
11 std::vector<int> hcsrRowPtrA = {0, 1, 4, 6, 8};
12 std::vector<int> hcsrColIndA = {0, 0, 1, 3, 1, 2, 0, 3};
13 std::vector<float> hcsrValA = {1.0f, 4.0f, 2.0f, 4.0f, 3.0f, 7.0f, 9.0f, 1.0f};
14
15 int* dcsrRowPtrA;
16 int* dcsrColIndA;
17 float* dcsrValA;
18 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrA, sizeof(int) * (m + 1)));
19 HIP_CHECK(hipMalloc((void**)&dcsrColIndA, sizeof(int) * nnz));
20 HIP_CHECK(hipMalloc((void**)&dcsrValA, sizeof(float) * nnz));
21
22 HIP_CHECK(
23 hipMemcpy(dcsrRowPtrA, hcsrRowPtrA.data(), sizeof(int) * (m + 1), hipMemcpyHostToDevice));
24 HIP_CHECK(hipMemcpy(dcsrColIndA, hcsrColIndA.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
25 HIP_CHECK(hipMemcpy(dcsrValA, hcsrValA.data(), sizeof(float) * nnz, hipMemcpyHostToDevice));
26
27 float* ddenseB;
28 HIP_CHECK(hipMalloc((void**)&ddenseB, sizeof(float) * m * n));
29
30 hipsparseHandle_t handle;
31 hipsparseSpMatDescr_t matA;
32 hipsparseDnMatDescr_t matB;
33
34 HIPSPARSE_CHECK(hipsparseCreate(&handle));
35
36 hipsparseIndexType_t rowIdxTypeA = HIPSPARSE_INDEX_32I;
37 hipsparseIndexType_t colIdxTypeA = HIPSPARSE_INDEX_32I;
38 hipDataType dataTypeA = HIP_R_32F;
39 hipsparseIndexBase_t idxBaseA = HIPSPARSE_INDEX_BASE_ZERO;
40
41 // Create sparse matrix A
42 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
43 m,
44 n,
45 nnz,
46 dcsrRowPtrA,
47 dcsrColIndA,
48 dcsrValA,
49 rowIdxTypeA,
50 colIdxTypeA,
51 idxBaseA,
52 dataTypeA));
53
54 // Create dense matrix B
55 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matB, m, n, m, ddenseB, HIP_R_32F, HIPSPARSE_ORDER_COL));
56
57 hipsparseSparseToDenseAlg_t alg = HIPSPARSE_SPARSETODENSE_ALG_DEFAULT;
58
59 size_t bufferSize;
60 HIPSPARSE_CHECK(hipsparseSparseToDense_bufferSize(handle, matA, matB, alg, &bufferSize));
61
62 void* tempBuffer;
63 HIP_CHECK(hipMalloc((void**)&tempBuffer, bufferSize));
64
65 // Complete the conversion
66 HIPSPARSE_CHECK(hipsparseSparseToDense(handle, matA, matB, alg, tempBuffer));
67
68 // Copy result back to host
69 std::vector<float> hdenseB(m * n);
70 HIP_CHECK(hipMemcpy(hdenseB.data(), ddenseB, sizeof(float) * m * n, hipMemcpyDeviceToHost));
71
72 // Clear hipSPARSE
73 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matA));
74 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matB));
75 HIPSPARSE_CHECK(hipsparseDestroy(handle));
76
77 // Clear device memory
78 HIP_CHECK(hipFree(dcsrRowPtrA));
79 HIP_CHECK(hipFree(dcsrColIndA));
80 HIP_CHECK(hipFree(dcsrValA));
81 HIP_CHECK(hipFree(ddenseB));
82 HIP_CHECK(hipFree(tempBuffer));
83
84 return 0;
85}
1int main(int argc, char* argv[])
2{
3 // 1 0 0 0
4 // A = 4 2 0 4
5 // 0 3 7 0
6 // 9 0 0 1
7 int m = 4;
8 int n = 4;
9 int nnz = 8;
10
11 int hcsrRowPtrA[] = {0, 1, 4, 6, 8};
12 int hcsrColIndA[] = {0, 0, 1, 3, 1, 2, 0, 3};
13 float hcsrValA[] = {1.0, 4.0, 2.0, 4.0, 3.0, 7.0, 9.0, 1.0};
14
15 int* dcsrRowPtrA;
16 int* dcsrColIndA;
17 float* dcsrValA;
18 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrA, sizeof(int) * (m + 1)));
19 HIP_CHECK(hipMalloc((void**)&dcsrColIndA, sizeof(int) * nnz));
20 HIP_CHECK(hipMalloc((void**)&dcsrValA, sizeof(float) * nnz));
21
22 HIP_CHECK(hipMemcpy(dcsrRowPtrA, hcsrRowPtrA, sizeof(int) * (m + 1), hipMemcpyHostToDevice));
23 HIP_CHECK(hipMemcpy(dcsrColIndA, hcsrColIndA, sizeof(int) * nnz, hipMemcpyHostToDevice));
24 HIP_CHECK(hipMemcpy(dcsrValA, hcsrValA, sizeof(float) * nnz, hipMemcpyHostToDevice));
25
26 float* ddenseB;
27 HIP_CHECK(hipMalloc((void**)&ddenseB, sizeof(float) * m * n));
28
29 hipsparseHandle_t handle;
30 hipsparseSpMatDescr_t matA;
31 hipsparseDnMatDescr_t matB;
32
33 HIPSPARSE_CHECK(hipsparseCreate(&handle));
34
35 hipsparseIndexType_t rowIdxTypeA = HIPSPARSE_INDEX_32I;
36 hipsparseIndexType_t colIdxTypeA = HIPSPARSE_INDEX_32I;
37 hipDataType dataTypeA = HIP_R_32F;
38 hipsparseIndexBase_t idxBaseA = HIPSPARSE_INDEX_BASE_ZERO;
39
40 // Create sparse matrix A
41 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
42 m,
43 n,
44 nnz,
45 dcsrRowPtrA,
46 dcsrColIndA,
47 dcsrValA,
48 rowIdxTypeA,
49 colIdxTypeA,
50 idxBaseA,
51 dataTypeA));
52
53 // Create dense matrix B
54 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matB, m, n, m, ddenseB, HIP_R_32F, HIPSPARSE_ORDER_COL));
55
56 hipsparseSparseToDenseAlg_t alg = HIPSPARSE_SPARSETODENSE_ALG_DEFAULT;
57
58 size_t bufferSize;
59 HIPSPARSE_CHECK(hipsparseSparseToDense_bufferSize(handle, matA, matB, alg, &bufferSize));
60
61 void* tempBuffer;
62 HIP_CHECK(hipMalloc((void**)&tempBuffer, bufferSize));
63
64 // Complete the conversion
65 HIPSPARSE_CHECK(hipsparseSparseToDense(handle, matA, matB, alg, tempBuffer));
66
67 // Copy result back to host
68 float* hdenseB = (float*)malloc((m * n) * sizeof(float));
69 HIP_CHECK(hipMemcpy(hdenseB, ddenseB, sizeof(float) * m * n, hipMemcpyDeviceToHost));
70
71 // Clear hipSPARSE
72 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matA));
73 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matB));
74 HIPSPARSE_CHECK(hipsparseDestroy(handle));
75
76 // Clear device memory
77 HIP_CHECK(hipFree(dcsrRowPtrA));
78 HIP_CHECK(hipFree(dcsrColIndA));
79 HIP_CHECK(hipFree(dcsrValA));
80 HIP_CHECK(hipFree(ddenseB));
81 HIP_CHECK(hipFree(tempBuffer));
82
83 return 0;
84}
hipsparseDenseToSparse_bufferSize()#
-
hipsparseStatus_t hipsparseDenseToSparse_bufferSize(hipsparseHandle_t handle, hipsparseConstDnMatDescr_t matA, hipsparseSpMatDescr_t matB, hipsparseDenseToSparseAlg_t alg, size_t *pBufferSizeInBytes)#
hipsparseDenseToSparse_bufferSizecomputes the required user allocated buffer size needed when converting a dense matrix to a sparse matrix. This routine currently accepts the sparse matrix descriptormatBin CSR, CSC, or COO format. This routine is used to determine the size of the buffer needed in hipsparseDenseToSparse_analysis and hipsparseDenseToSparse_convert.hipsparseDenseToSparse_bufferSizesupports different data types for the dense and sparse matrices. See hipsparseDenseToSparse_convert for a complete listing of all the data types available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
matA – [in] dense matrix descriptor.
matB – [in] sparse matrix descriptor.
alg – [in] algorithm for the dense to sparse computation.
pBufferSizeInBytes – [out] number of bytes of the temporary storage buffer.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,matA,matB, orpBufferSizeInBytespointer is invalid.
hipsparseDenseToSparse_analysis()#
-
hipsparseStatus_t hipsparseDenseToSparse_analysis(hipsparseHandle_t handle, hipsparseConstDnMatDescr_t matA, hipsparseSpMatDescr_t matB, hipsparseDenseToSparseAlg_t alg, void *externalBuffer)#
hipsparseDenseToSparse_analysisperforms analysis that is later used in hipsparseDenseToSparse_convert when converting a dense matrix to sparse matrix. This routine currently accepts the sparse matrix descriptormatBin CSR, CSC, or COO format. This routine takes a user allocated buffer whose size must first be computed using hipsparseDenseToSparse_bufferSize.hipsparseDenseToSparse_analysissupports different data types for the dense and sparse matrices. See hipsparseDenseToSparse_convert for a complete listing of all the data types available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
matA – [in] dense matrix descriptor.
matB – [in] sparse matrix descriptor.
alg – [in] algorithm for the dense to sparse computation.
externalBuffer – [out] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,matA,matB, orexternalBufferpointer is invalid.
hipsparseDenseToSparse_convert()#
-
hipsparseStatus_t hipsparseDenseToSparse_convert(hipsparseHandle_t handle, hipsparseConstDnMatDescr_t matA, hipsparseSpMatDescr_t matB, hipsparseDenseToSparseAlg_t alg, void *externalBuffer)#
Dense matrix to sparse matrix conversion.
hipsparseDenseToSparse_convertconverts a dense matrix to a sparse matrix. This routine currently accepts the sparse matrix descriptormatBin CSR, CSC, or COO format. This routine requires a user allocated buffer whose size must be determined by first calling hipsparseDenseToSparse_bufferSize.The conversion of a dense matrix into a sparse one involves three steps. First, the user creates the dense and sparse matrix descriptors. Because the number of non-zeros that will exist in the sparse matrix is not known apriori, when creating the sparse matrix descriptor, the user simply sets the arrays to
NULLand the non-zero count to zero. For example, in the case of a CSR sparse matrix, this would look like:In the case of a COO sparse matrix, this would look like:hipsparseCreateCsr(&matB, m, n, 0, dcsrRowPtrB, // This array can be allocated as its size (i.e. m + 1) is known NULL, // Column indices array size is not yet known, pass NULL for now NULL, // Values array size is not yet known, pass NULL for now rowIdxTypeB, colIdxTypeB, idxBaseB, dataTypeB);
Once the descriptors have been created, the user calls hipsparseDenseToSparse_bufferSize. This routine will determine the size of the required temporary storage buffer. The user then allocates this buffer and passes it to hipsparseDenseToSparse_analysis which will perform analysis on the dense matrix in order to determine the number of non-zeros that will exist in the sparse matrix. Once this hipsparseDenseToSparse_analysis routine has been called, the non-zero count is stored in the sparse matrix descriptorhipsparseCreateCoo(&matB, m, n, 0, NULL, // Row indices array size is not yet known, pass NULL for now NULL, // Column indices array size is not yet known, pass NULL for now NULL, // Values array size is not yet known, pass NULL for now rowIdxTypeB, colIdxTypeB, idxBaseB, dataTypeB);
matB. In order to allocate our remaining sparse matrix arrays, we query the sparse matrix descriptormatBfor this non-zero count:The remaining arrays are then allocated and set on the sparse matrix descriptor// Grab the non-zero count from the B matrix decriptor int64_t rows; int64_t cols; int64_t nnz; hipsparseSpMatGetSize(matB, &rows, &cols, &nnz);
matB. Finally, we complete the conversion by calling hipsparseDenseToSparse_convert. Once the conversion is complete, the user is free to deallocate the storage buffer. See full example below for details.hipsparseDenseToSparse_convertsupports the following uniform precision data types for the dense and sparse matrices \(A\) and \(B\):- Uniform Precisions:
A / B
HIP_R_16F
HIP_R_16BF
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
Note
Currently only the sparse matrix formats CSR, CSC, and COO are supported when converting a dense matrix to a sparse matrix.
- Parameters:
handle – [in] handle to the hipsparse library context queue.
matA – [in] dense matrix descriptor.
matB – [in] sparse matrix descriptor.
alg – [in] algorithm for the dense to sparse computation.
externalBuffer – [out] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,matA,matB, orexternalBufferpointer is invalid.
1int main(int argc, char* argv[])
2{
3 // 1 0 0 0
4 // A = 4 2 0 4
5 // 0 3 7 0
6 // 9 0 0 1
7 int m = 4;
8 int n = 4;
9
10 std::vector<float> hdenseA = {1.0f,
11 4.0f,
12 0.0f,
13 9.0f,
14 0.0f,
15 2.0f,
16 3.0f,
17 0.0f,
18 0.0f,
19 0.0f,
20 7.0f,
21 0.0f,
22 0.0f,
23 4.0f,
24 0.0f,
25 1.0f};
26
27 float* ddenseA;
28 HIP_CHECK(hipMalloc((void**)&ddenseA, sizeof(float) * m * n));
29 HIP_CHECK(hipMemcpy(ddenseA, hdenseA.data(), sizeof(float) * m * n, hipMemcpyHostToDevice));
30
31 int* dcsrRowPtrB;
32 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrB, sizeof(int) * (m + 1)));
33
34 hipsparseHandle_t handle;
35 hipsparseDnMatDescr_t matA;
36 hipsparseSpMatDescr_t matB;
37
38 HIPSPARSE_CHECK(hipsparseCreate(&handle));
39
40 // Create dense matrix A
41 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matA, m, n, m, ddenseA, HIP_R_32F, HIPSPARSE_ORDER_COL));
42
43 hipsparseIndexType_t rowIdxTypeB = HIPSPARSE_INDEX_32I;
44 hipsparseIndexType_t colIdxTypeB = HIPSPARSE_INDEX_32I;
45 hipDataType dataTypeB = HIP_R_32F;
46 hipsparseIndexBase_t idxBaseB = HIPSPARSE_INDEX_BASE_ZERO;
47
48 // Create sparse matrix B
49 HIPSPARSE_CHECK(hipsparseCreateCsr(
50 &matB, m, n, 0, dcsrRowPtrB, NULL, NULL, rowIdxTypeB, colIdxTypeB, idxBaseB, dataTypeB));
51
52 hipsparseDenseToSparseAlg_t alg = HIPSPARSE_DENSETOSPARSE_ALG_DEFAULT;
53
54 size_t bufferSize;
55 HIPSPARSE_CHECK(hipsparseDenseToSparse_bufferSize(handle, matA, matB, alg, &bufferSize));
56
57 void* tempBuffer;
58 HIP_CHECK(hipMalloc((void**)&tempBuffer, bufferSize));
59
60 // Perform analysis which will determine the number of non-zeros in the CSR matrix
61 HIPSPARSE_CHECK(hipsparseDenseToSparse_analysis(handle, matA, matB, alg, tempBuffer));
62
63 // Grab the non-zero count from the B matrix decriptor
64 int64_t rows;
65 int64_t cols;
66 int64_t nnz;
67 HIPSPARSE_CHECK(hipsparseSpMatGetSize(matB, &rows, &cols, &nnz));
68
69 // Allocate the column indices and values arrays
70 int* dcsrColIndB;
71 float* dcsrValB;
72 HIP_CHECK(hipMalloc((void**)&dcsrColIndB, sizeof(int) * nnz));
73 HIP_CHECK(hipMalloc((void**)&dcsrValB, sizeof(float) * nnz));
74
75 // Set the newly allocated arrays on the sparse matrix descriptor
76 HIPSPARSE_CHECK(hipsparseCsrSetPointers(matB, dcsrRowPtrB, dcsrColIndB, dcsrValB));
77
78 // Complete the conversion
79 HIPSPARSE_CHECK(hipsparseDenseToSparse_convert(handle, matA, matB, alg, tempBuffer));
80
81 // Copy result back to host
82 std::vector<int> hcsrRowPtrB(m + 1);
83 std::vector<int> hcsrColIndB(nnz);
84 std::vector<float> hcsrValB(nnz);
85 HIP_CHECK(
86 hipMemcpy(hcsrRowPtrB.data(), dcsrRowPtrB, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
87 HIP_CHECK(hipMemcpy(hcsrColIndB.data(), dcsrColIndB, sizeof(int) * nnz, hipMemcpyDeviceToHost));
88 HIP_CHECK(hipMemcpy(hcsrValB.data(), dcsrValB, sizeof(float) * nnz, hipMemcpyDeviceToHost));
89
90 // Clear hipSPARSE
91 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matA));
92 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matB));
93 HIPSPARSE_CHECK(hipsparseDestroy(handle));
94
95 // Clear device memory
96 HIP_CHECK(hipFree(ddenseA));
97 HIP_CHECK(hipFree(dcsrRowPtrB));
98 HIP_CHECK(hipFree(dcsrColIndB));
99 HIP_CHECK(hipFree(dcsrValB));
100 HIP_CHECK(hipFree(tempBuffer));
101
102 return 0;
103}
1int main(int argc, char* argv[])
2{
3 // 1 0 0 0
4 // A = 4 2 0 4
5 // 0 3 7 0
6 // 9 0 0 1
7 int m = 4;
8 int n = 4;
9
10 float hdenseA[]
11 = {1.0, 4.0, 0.0, 9.0, 0.0, 2.0, 3.0, 0.0, 0.0, 0.0, 7.0, 0.0, 0.0, 4.0, 0.0, 1.0};
12
13 float* ddenseA;
14 HIP_CHECK(hipMalloc((void**)&ddenseA, sizeof(float) * m * n));
15 HIP_CHECK(hipMemcpy(ddenseA, hdenseA, sizeof(float) * m * n, hipMemcpyHostToDevice));
16
17 int* dcsrRowPtrB;
18 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrB, sizeof(int) * (m + 1)));
19
20 hipsparseHandle_t handle;
21 hipsparseDnMatDescr_t matA;
22 hipsparseSpMatDescr_t matB;
23
24 HIPSPARSE_CHECK(hipsparseCreate(&handle));
25
26 // Create dense matrix A
27 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matA, m, n, m, ddenseA, HIP_R_32F, HIPSPARSE_ORDER_COL));
28
29 hipsparseIndexType_t rowIdxTypeB = HIPSPARSE_INDEX_32I;
30 hipsparseIndexType_t colIdxTypeB = HIPSPARSE_INDEX_32I;
31 hipDataType dataTypeB = HIP_R_32F;
32 hipsparseIndexBase_t idxBaseB = HIPSPARSE_INDEX_BASE_ZERO;
33
34 // Create sparse matrix B
35 HIPSPARSE_CHECK(hipsparseCreateCsr(
36 &matB, m, n, 0, dcsrRowPtrB, NULL, NULL, rowIdxTypeB, colIdxTypeB, idxBaseB, dataTypeB));
37
38 hipsparseDenseToSparseAlg_t alg = HIPSPARSE_DENSETOSPARSE_ALG_DEFAULT;
39
40 size_t bufferSize;
41 HIPSPARSE_CHECK(hipsparseDenseToSparse_bufferSize(handle, matA, matB, alg, &bufferSize));
42
43 void* tempBuffer;
44 HIP_CHECK(hipMalloc((void**)&tempBuffer, bufferSize));
45
46 // Perform analysis which will determine the number of non-zeros in the CSR matrix
47 HIPSPARSE_CHECK(hipsparseDenseToSparse_analysis(handle, matA, matB, alg, tempBuffer));
48
49 // Grab the non-zero count from the B matrix decriptor
50 int64_t rows;
51 int64_t cols;
52 int64_t nnz;
53 HIPSPARSE_CHECK(hipsparseSpMatGetSize(matB, &rows, &cols, &nnz));
54
55 // Allocate the column indices and values arrays
56 int* dcsrColIndB;
57 float* dcsrValB;
58 HIP_CHECK(hipMalloc((void**)&dcsrColIndB, sizeof(int) * nnz));
59 HIP_CHECK(hipMalloc((void**)&dcsrValB, sizeof(float) * nnz));
60
61 // Set the newly allocated arrays on the sparse matrix descriptor
62 HIPSPARSE_CHECK(hipsparseCsrSetPointers(matB, dcsrRowPtrB, dcsrColIndB, dcsrValB));
63
64 // Complete the conversion
65 HIPSPARSE_CHECK(hipsparseDenseToSparse_convert(handle, matA, matB, alg, tempBuffer));
66
67 // Copy result back to host
68 int* hcsrRowPtrB = (int*)malloc((m + 1) * sizeof(int));
69 int* hcsrColIndB = (int*)malloc((nnz) * sizeof(int));
70 float hcsrValB[nnz];
71 HIP_CHECK(hipMemcpy(hcsrRowPtrB, dcsrRowPtrB, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
72 HIP_CHECK(hipMemcpy(hcsrColIndB, dcsrColIndB, sizeof(int) * nnz, hipMemcpyDeviceToHost));
73 HIP_CHECK(hipMemcpy(hcsrValB, dcsrValB, sizeof(float) * nnz, hipMemcpyDeviceToHost));
74
75 // Clear hipSPARSE
76 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matA));
77 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matB));
78 HIPSPARSE_CHECK(hipsparseDestroy(handle));
79
80 // Clear device memory
81 HIP_CHECK(hipFree(ddenseA));
82 HIP_CHECK(hipFree(dcsrRowPtrB));
83 HIP_CHECK(hipFree(dcsrColIndB));
84 HIP_CHECK(hipFree(dcsrValB));
85 HIP_CHECK(hipFree(tempBuffer));
86
87 return 0;
88}
hipsparseSpVV_bufferSize()#
-
hipsparseStatus_t hipsparseSpVV_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t opX, hipsparseConstSpVecDescr_t vecX, hipsparseConstDnVecDescr_t vecY, void *result, hipDataType computeType, size_t *pBufferSizeInBytes)#
hipsparseSpVV_bufferSizecomputes the required user allocated buffer size needed when computing the inner dot product of a sparse vector with a dense vector:\[ \text{result} := op(x) \cdot y, \]hipsparseSpVV_bufferSizesupports multiple combinations of data types and compute types. See hipsparseSpVV for a complete listing of all the data type and compute type combinations available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opX – [in] sparse vector operation type.
vecX – [in] sparse vector descriptor.
vecY – [in] dense vector descriptor.
result – [out] pointer to the result, can be host or device memory.
computeType – [in] floating point precision for the SpVV computation.
pBufferSizeInBytes – [out] number of bytes of the temporary storage buffer.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_NOT_INITIALIZED –
handleis not initialized.HIPSPARSE_STATUS_INVALID_VALUE –
handle,vecX,vecY,resultorpBufferSizeInBytesis nullptr.HIPSPARSE_STATUS_NOT_SUPPORTED –
computeTypeis currently not supported.
hipsparseSpVV()#
-
hipsparseStatus_t hipsparseSpVV(hipsparseHandle_t handle, hipsparseOperation_t opX, hipsparseConstSpVecDescr_t vecX, hipsparseConstDnVecDescr_t vecY, void *result, hipDataType computeType, void *externalBuffer)#
Compute the inner dot product of a sparse vector with a dense vector.
hipsparseSpVVcomputes the inner dot product of the sparse vector \(x\) with the dense vector \(y\), such that\[ \text{result} := op(x) \cdot y, \]with\[\begin{split} op(x) = \left\{ \begin{array}{ll} x, & \text{if trans == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ \bar{x}, & \text{if trans == HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE} \\ \end{array} \right. \end{split}\]result = 0; for(i = 0; i < nnz; ++i) { result += x_val[i] * y[x_ind[i]]; }
Performing the above operation involves two steps. First, the user calls
hipsparseSpVV_bufferSizewhich will return the required temporary buffer size. The user then allocates this buffer. Finally, the user then completes the computation by callinghipsparseSpVVwith the newly allocated buffer. Once the computation is complete, the user is free to deallocate the buffer.hipsparseSpVVsupports the following uniform and mixed precision data types for the sparse and dense vectors \(x\) and \(y\) and compute types for the scalar \(result\).- Uniform Precisions:
X / Y / compute_type
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
- Mixed precisions:
X / Y
compute_type
HIP_R_8I
HIP_R_32I
HIP_R_8I
HIP_R_32F
HIP_R_16F
HIP_R_32F
HIP_R_16BF
HIP_R_32F
- Parameters:
handle – [in] handle to the hipsparse library context queue.
opX – [in] sparse vector operation type.
vecX – [in] sparse vector descriptor.
vecY – [in] dense vector descriptor.
result – [out] pointer to the result, can be host or device memory
computeType – [in] floating point precision for the SpVV computation.
externalBuffer – [out] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,vecX,vecY,resultorexternalBufferpointer is invalid.HIPSPARSE_STATUS_NOT_SUPPORTED –
computeTypeis currently not supported.
1int main(int argc, char* argv[])
2{
3 // Number of non-zeros of the sparse vector
4 int nnz = 3;
5
6 // Size of sparse and dense vector
7 int size = 9;
8
9 // Sparse index vector
10 std::vector<int> hxInd = {0, 3, 5};
11
12 // Sparse value vector
13 std::vector<float> hxVal = {1.0f, 2.0f, 3.0f};
14
15 // Dense vector
16 std::vector<float> hy = {1.0f, 2.0f, 3.0f, 4.0f, 5.0f, 6.0f, 7.0f, 8.0f, 9.0f};
17
18 // Offload data to device
19 int* dxInd;
20 float* dxVal;
21 float* dy;
22 HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
23 HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
24 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
25
26 HIP_CHECK(hipMemcpy(dxInd, hxInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
27 HIP_CHECK(hipMemcpy(dxVal, hxVal.data(), sizeof(float) * nnz, hipMemcpyHostToDevice));
28 HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(float) * size, hipMemcpyHostToDevice));
29
30 hipsparseHandle_t handle;
31 HIPSPARSE_CHECK(hipsparseCreate(&handle));
32
33 // Create sparse vector X
34 hipsparseSpVecDescr_t vecX;
35 HIPSPARSE_CHECK(hipsparseCreateSpVec(
36 &vecX, size, nnz, dxInd, dxVal, HIPSPARSE_INDEX_32I, HIPSPARSE_INDEX_BASE_ZERO, HIP_R_32F));
37
38 // Create dense vector Y
39 hipsparseDnVecDescr_t vecY;
40 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, size, dy, HIP_R_32F));
41
42 // Obtain buffer size
43 float hresult = 0.0f;
44 size_t buffer_size;
45 HIPSPARSE_CHECK(hipsparseSpVV_bufferSize(
46 handle, HIPSPARSE_OPERATION_NON_TRANSPOSE, vecX, vecY, &hresult, HIP_R_32F, &buffer_size));
47
48 void* temp_buffer;
49 HIP_CHECK(hipMalloc(&temp_buffer, buffer_size));
50
51 // SpVV
52 HIPSPARSE_CHECK(hipsparseSpVV(
53 handle, HIPSPARSE_OPERATION_NON_TRANSPOSE, vecX, vecY, &hresult, HIP_R_32F, temp_buffer));
54
55 HIP_CHECK(hipDeviceSynchronize());
56
57 std::cout << "hresult: " << hresult << std::endl;
58
59 // Clear hipSPARSE
60 HIPSPARSE_CHECK(hipsparseDestroySpVec(vecX));
61 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
62 HIPSPARSE_CHECK(hipsparseDestroy(handle));
63
64 // Clear device memory
65 HIP_CHECK(hipFree(dxInd));
66 HIP_CHECK(hipFree(dxVal));
67 HIP_CHECK(hipFree(dy));
68 HIP_CHECK(hipFree(temp_buffer));
69
70 return 0;
71}
1int main(int argc, char* argv[])
2{
3 // Number of non-zeros of the sparse vector
4 int nnz = 3;
5
6 // Size of sparse and dense vector
7 int size = 9;
8
9 // Sparse index vector
10 int hxInd[] = {0, 3, 5};
11
12 // Sparse value vector
13 float hxVal[] = {1.0, 2.0, 3.0};
14
15 // Dense vector
16 float hy[] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0};
17
18 // Offload data to device
19 int* dxInd;
20 float* dxVal;
21 float* dy;
22 HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
23 HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
24 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
25
26 HIP_CHECK(hipMemcpy(dxInd, hxInd, sizeof(int) * nnz, hipMemcpyHostToDevice));
27 HIP_CHECK(hipMemcpy(dxVal, hxVal, sizeof(float) * nnz, hipMemcpyHostToDevice));
28 HIP_CHECK(hipMemcpy(dy, hy, sizeof(float) * size, hipMemcpyHostToDevice));
29
30 hipsparseHandle_t handle;
31 HIPSPARSE_CHECK(hipsparseCreate(&handle));
32
33 // Create sparse vector X
34 hipsparseSpVecDescr_t vecX;
35 HIPSPARSE_CHECK(hipsparseCreateSpVec(
36 &vecX, size, nnz, dxInd, dxVal, HIPSPARSE_INDEX_32I, HIPSPARSE_INDEX_BASE_ZERO, HIP_R_32F));
37
38 // Create dense vector Y
39 hipsparseDnVecDescr_t vecY;
40 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, size, dy, HIP_R_32F));
41
42 // Obtain buffer size
43 float hresult = 0.0;
44 size_t buffer_size;
45 HIPSPARSE_CHECK(hipsparseSpVV_bufferSize(
46 handle, HIPSPARSE_OPERATION_NON_TRANSPOSE, vecX, vecY, &hresult, HIP_R_32F, &buffer_size));
47
48 void* temp_buffer;
49 HIP_CHECK(hipMalloc(&temp_buffer, buffer_size));
50
51 // SpVV
52 HIPSPARSE_CHECK(hipsparseSpVV(
53 handle, HIPSPARSE_OPERATION_NON_TRANSPOSE, vecX, vecY, &hresult, HIP_R_32F, temp_buffer));
54
55 HIP_CHECK(hipDeviceSynchronize());
56
57 printf("hresult: %f\n", hresult);
58
59 // Clear hipSPARSE
60 HIPSPARSE_CHECK(hipsparseDestroySpVec(vecX));
61 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
62 HIPSPARSE_CHECK(hipsparseDestroy(handle));
63
64 // Clear device memory
65 HIP_CHECK(hipFree(dxInd));
66 HIP_CHECK(hipFree(dxVal));
67 HIP_CHECK(hipFree(dy));
68 HIP_CHECK(hipFree(temp_buffer));
69
70 return 0;
71}
hipsparseSpMV_bufferSize()#
-
hipsparseStatus_t hipsparseSpMV_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t opA, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstDnVecDescr_t vecX, const void *beta, const hipsparseDnVecDescr_t vecY, hipDataType computeType, hipsparseSpMVAlg_t alg, size_t *pBufferSizeInBytes)#
hipsparseSpMV_bufferSizecomputes the required user allocated buffer size needed when computing the sparse matrix multiplication with a dense vector:\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]where \(op(A)\) is a sparse \(m \times n\) matrix in CSR, CSC, COO, or COO (AoS) format, \(x\) is a dense vector of length \(n\) and \(y\) is a dense vector of length \(m\).hipsparseSpMV_bufferSizesupports multiple combinations of data types and compute types. See hipsparseSpMV for a complete listing of all the data type and compute type combinations available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] matrix operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] sparse matrix descriptor.
vecX – [in] dense vector descriptor.
beta – [in] scalar \(\beta\).
vecY – [inout] dense vector descriptor.
computeType – [in] floating point precision for the SpMV computation.
alg – [in] SpMV algorithm for the SpMV computation.
pBufferSizeInBytes – [out] number of bytes of the temporary storage buffer.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_NOT_INITIALIZED –
handleis not initialized.HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,matA,vecX,beta,vecYorpBufferSizeInBytesis nullptr,opAis invalid, or the matrix or vector dimensions are incompatible.HIPSPARSE_STATUS_NOT_SUPPORTED –
computeTypeoralgis currently not supported.
hipsparseSpMV_preprocess()#
-
hipsparseStatus_t hipsparseSpMV_preprocess(hipsparseHandle_t handle, hipsparseOperation_t opA, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstDnVecDescr_t vecX, const void *beta, const hipsparseDnVecDescr_t vecY, hipDataType computeType, hipsparseSpMVAlg_t alg, void *externalBuffer)#
hipsparseSpMV_preprocessperforms analysis on the sparse matrix \(op(A)\) when computing the sparse matrix multiplication with a dense vector:\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]where \(op(A)\) is a sparse \(m \times n\) matrix in CSR, CSC, COO, or COO (AoS) format, \(x\) is a dense vector of length \(n\) and \(y\) is a dense vector of length \(m\). This step is optional but if used may results in better performance.hipsparseSpMV_preprocesssupports multiple combinations of data types and compute types. See hipsparseSpMV for a complete listing of all the data type and compute type combinations available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] matrix operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] matrix descriptor.
vecX – [in] vector descriptor.
beta – [in] scalar \(\beta\).
vecY – [inout] vector descriptor.
computeType – [in] floating point precision for the SpMV computation.
alg – [in] SpMV algorithm for the SpMV computation.
externalBuffer – [out] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,matA,x,beta,yorexternalBufferpointer is invalid or ifopA,computeType,algis incorrect.HIPSPARSE_STATUS_NOT_SUPPORTED –
computeTypeoralgis currently not supported.
hipsparseSpMV()#
-
hipsparseStatus_t hipsparseSpMV(hipsparseHandle_t handle, hipsparseOperation_t opA, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstDnVecDescr_t vecX, const void *beta, const hipsparseDnVecDescr_t vecY, hipDataType computeType, hipsparseSpMVAlg_t alg, void *externalBuffer)#
Compute the sparse matrix multiplication with a dense vector.
hipsparseSpMVmultiplies the scalar \(\alpha\) with a sparse \(m \times n\) matrix \(op(A)\), defined in CSR, CSC, COO, or COO (AoS) format, with the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == 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 involves multiple steps. First the user calls hipsparseSpMV_bufferSize to determine the size of the required temporary storage buffer. The user then allocates this buffer and calls hipsparseSpMV_preprocess. Depending on the algorithm and sparse matrix format, this will perform analysis on the sparsity pattern of \(op(A)\). Finally the user completes the operation by calling
hipsparseSpMV. The buffer size and preprecess routines only need to be called once for a given sparse matrix \(op(A)\) while the computation can be repeatedly used with different \(x\) and \(y\) vectors. Once all calls tohipsparseSpMVare complete, the temporary buffer can be deallocated.hipsparseSpMVsupports multiple different algorithms. These algorithms have different trade offs depending on the sparsity pattern of the matrix, whether or not the results need to be deterministic, and how many times the sparse-vector product will be performed.CSR Algorithms
HIPSPARSE_SPMV_CSR_ALG1
HIPSPARSE_SPMV_CSR_ALG2
COO Algorithms
HIPSPARSE_SPMV_COO_ALG1
HIPSPARSE_SPMV_COO_ALG2
hipsparseSpMVsupports multiple combinations of data types and compute types. The tables below indicate the currently supported data types that can be used for the sparse matrix \(op(A)\) and the dense vectors \(x\) and \(y\) and the compute type for \(\alpha\) and \(\beta\). The advantage of using different data types is to save on memory bandwidth and storage when a user application allows while performing the actual computation in a higher precision.hipsparseSpMVsupports HIPSPARSE_INDEX_32I and HIPSPARSE_INDEX_64I index precisions for storing the row pointer and row/column indices arrays of the sparse matrices.- Uniform Precisions:
A / X / Y / compute_type
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
- Mixed precisions:
A / X
Y
compute_type
HIP_R_8I
HIP_R_32I
HIP_R_32I
HIP_R_8I
HIP_R_32F
HIP_R_32F
HIP_R_16F
HIP_R_32F
HIP_R_32F
HIP_R_16BF
HIP_R_32F
HIP_R_32F
- Mixed-regular real precisions
A
X / Y / compute_type
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
- Mixed-regular Complex precisions
A
X / Y / compute_type
HIP_R_32F
HIP_C_32F
HIP_R_64F
HIP_C_64F
Note
None of the algorithms above are deterministic when \(A\) is transposed.
Note
The sparse matrix formats currently supported are: HIPSPARSE_FORMAT_COO, HIPSPARSE_FORMAT_COO_AOS, HIPSPARSE_FORMAT_CSR, and HIPSPARSE_FORMAT_CSC.
Note
Only the hipsparseSpMV_bufferSize and hipsparseSpMV routines are non blocking and executed asynchronously with respect to the host. They may return before the actual computation has finished. The hipsparseSpMV_preprocess routine is blocking with respect to the host.
Note
Only the hipsparseSpMV_bufferSize and the hipsparseSpMV routines support execution in a hipGraph context. The hipsparseSpMV_preprocess stage does not support hipGraph.
- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] matrix operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] matrix descriptor.
vecX – [in] vector descriptor.
beta – [in] scalar \(\beta\).
vecY – [inout] vector descriptor.
computeType – [in] floating point precision for the SpMV computation.
alg – [in] SpMV algorithm for the SpMV computation.
externalBuffer – [out] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,matA,x,beta,yorexternalBufferpointer is invalid or ifopA,computeType,algis incorrect.HIPSPARSE_STATUS_NOT_SUPPORTED –
computeTypeoralgis currently not supported.
1int main(int argc, char* argv[])
2{
3 // A, x, and y are m×k, k×1, and m×1
4 int m = 3, k = 4;
5 int nnz_A = 8;
6 hipsparseOperation_t transA = HIPSPARSE_OPERATION_NON_TRANSPOSE;
7
8 // alpha and beta
9 float alpha = 0.5f;
10 float beta = 0.25f;
11
12 std::vector<int> hcsrRowPtr = {0, 3, 5, 8};
13 std::vector<int> hcsrColInd = {0, 1, 3, 1, 2, 0, 2, 3};
14 std::vector<float> hcsrVal = {1, 2, 3, 4, 5, 6, 7, 8};
15
16 std::vector<float> hx(k, 1.0f);
17 std::vector<float> hy(m, 1.0f);
18
19 int* dcsrRowPtr;
20 int* dcsrColInd;
21 float* dcsrVal;
22 HIP_CHECK(hipMalloc((void**)&dcsrRowPtr, sizeof(int) * (m + 1)));
23 HIP_CHECK(hipMalloc((void**)&dcsrColInd, sizeof(int) * nnz_A));
24 HIP_CHECK(hipMalloc((void**)&dcsrVal, sizeof(float) * nnz_A));
25
26 HIP_CHECK(
27 hipMemcpy(dcsrRowPtr, hcsrRowPtr.data(), sizeof(int) * (m + 1), hipMemcpyHostToDevice));
28 HIP_CHECK(hipMemcpy(dcsrColInd, hcsrColInd.data(), sizeof(int) * nnz_A, hipMemcpyHostToDevice));
29 HIP_CHECK(hipMemcpy(dcsrVal, hcsrVal.data(), sizeof(float) * nnz_A, hipMemcpyHostToDevice));
30
31 hipsparseHandle_t handle;
32 HIPSPARSE_CHECK(hipsparseCreate(&handle));
33
34 hipsparseSpMatDescr_t matA;
35 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
36 m,
37 k,
38 nnz_A,
39 dcsrRowPtr,
40 dcsrColInd,
41 dcsrVal,
42 HIPSPARSE_INDEX_32I,
43 HIPSPARSE_INDEX_32I,
44 HIPSPARSE_INDEX_BASE_ZERO,
45 HIP_R_32F));
46
47 // Allocate memory for the vector x
48 float* dx;
49 HIP_CHECK(hipMalloc((void**)&dx, sizeof(float) * k));
50 HIP_CHECK(hipMemcpy(dx, hx.data(), sizeof(float) * k, hipMemcpyHostToDevice));
51
52 hipsparseDnVecDescr_t vecX;
53 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecX, k, dx, HIP_R_32F));
54
55 // Allocate memory for the resulting vector y
56 float* dy;
57 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * m));
58 HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(float) * m, hipMemcpyHostToDevice));
59
60 hipsparseDnMatDescr_t vecY;
61 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, m, dy, HIP_R_32F));
62
63 // Compute buffersize
64 size_t bufferSize;
65 HIPSPARSE_CHECK(hipsparseSpMV_bufferSize(handle,
66 transA,
67 &alpha,
68 matA,
69 vecX,
70 &beta,
71 vecY,
72 HIP_R_32F,
73 HIPSPARSE_MV_ALG_DEFAULT,
74 &bufferSize));
75
76 void* buffer;
77 HIP_CHECK(hipMalloc(&buffer, bufferSize));
78
79 // Preprocess operation (Optional)
80 HIPSPARSE_CHECK(hipsparseSpMV_preprocess(handle,
81 transA,
82 &alpha,
83 matA,
84 vecX,
85 &beta,
86 vecY,
87 HIP_R_32F,
88 HIPSPARSE_MV_ALG_DEFAULT,
89 buffer));
90
91 // Perform operation
92 HIPSPARSE_CHECK(hipsparseSpMV(handle,
93 transA,
94 &alpha,
95 matA,
96 vecX,
97 &beta,
98 vecY,
99 HIP_R_32F,
100 HIPSPARSE_MV_ALG_DEFAULT,
101 buffer));
102
103 // Copy device to host
104 HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(float) * m, hipMemcpyDeviceToHost));
105
106 // Destroy matrix descriptors and handles
107 HIPSPARSE_CHECK(hipsparseDestroySpMat(matA));
108 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecX));
109 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
110 HIPSPARSE_CHECK(hipsparseDestroy(handle));
111
112 HIP_CHECK(hipFree(buffer));
113 HIP_CHECK(hipFree(dcsrRowPtr));
114 HIP_CHECK(hipFree(dcsrColInd));
115 HIP_CHECK(hipFree(dcsrVal));
116 HIP_CHECK(hipFree(dx));
117 HIP_CHECK(hipFree(dy));
118
119 return 0;
120}
1int main(int argc, char* argv[])
2{
3 // A, x, and y are m×k, k×1, and m×1
4 int m = 3, k = 4;
5 int nnz_A = 8;
6 hipsparseOperation_t transA = HIPSPARSE_OPERATION_NON_TRANSPOSE;
7
8 // alpha and beta
9 float alpha = 0.5;
10 float beta = 0.25;
11
12 int hcsrRowPtr[] = {0, 3, 5, 8};
13 int hcsrColInd[] = {0, 1, 3, 1, 2, 0, 2, 3};
14 float hcsrVal[] = {1, 2, 3, 4, 5, 6, 7, 8};
15
16 float* hx = (float*)malloc(k * sizeof(float));
17 float* hy = (float*)malloc(m * sizeof(float));
18
19 int* dcsrRowPtr;
20 int* dcsrColInd;
21 float* dcsrVal;
22 HIP_CHECK(hipMalloc((void**)&dcsrRowPtr, sizeof(int) * (m + 1)));
23 HIP_CHECK(hipMalloc((void**)&dcsrColInd, sizeof(int) * nnz_A));
24 HIP_CHECK(hipMalloc((void**)&dcsrVal, sizeof(float) * nnz_A));
25
26 HIP_CHECK(hipMemcpy(dcsrRowPtr, hcsrRowPtr, sizeof(int) * (m + 1), hipMemcpyHostToDevice));
27 HIP_CHECK(hipMemcpy(dcsrColInd, hcsrColInd, sizeof(int) * nnz_A, hipMemcpyHostToDevice));
28 HIP_CHECK(hipMemcpy(dcsrVal, hcsrVal, sizeof(float) * nnz_A, hipMemcpyHostToDevice));
29
30 hipsparseHandle_t handle;
31 HIPSPARSE_CHECK(hipsparseCreate(&handle));
32
33 hipsparseSpMatDescr_t matA;
34 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
35 m,
36 k,
37 nnz_A,
38 dcsrRowPtr,
39 dcsrColInd,
40 dcsrVal,
41 HIPSPARSE_INDEX_32I,
42 HIPSPARSE_INDEX_32I,
43 HIPSPARSE_INDEX_BASE_ZERO,
44 HIP_R_32F));
45
46 // Allocate memory for the vector x
47 float* dx;
48 HIP_CHECK(hipMalloc((void**)&dx, sizeof(float) * k));
49 HIP_CHECK(hipMemcpy(dx, hx, sizeof(float) * k, hipMemcpyHostToDevice));
50
51 hipsparseDnVecDescr_t vecX;
52 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecX, k, dx, HIP_R_32F));
53
54 // Allocate memory for the resulting vector y
55 float* dy;
56 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * m));
57 HIP_CHECK(hipMemcpy(dy, hy, sizeof(float) * m, hipMemcpyHostToDevice));
58
59 hipsparseDnMatDescr_t vecY;
60 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, m, dy, HIP_R_32F));
61
62 // Compute buffersize
63 size_t bufferSize;
64 HIPSPARSE_CHECK(hipsparseSpMV_bufferSize(handle,
65 transA,
66 &alpha,
67 matA,
68 vecX,
69 &beta,
70 vecY,
71 HIP_R_32F,
72 HIPSPARSE_MV_ALG_DEFAULT,
73 &bufferSize));
74
75 void* buffer;
76 HIP_CHECK(hipMalloc(&buffer, bufferSize));
77
78 // Preprocess operation (Optional)
79 HIPSPARSE_CHECK(hipsparseSpMV_preprocess(handle,
80 transA,
81 &alpha,
82 matA,
83 vecX,
84 &beta,
85 vecY,
86 HIP_R_32F,
87 HIPSPARSE_MV_ALG_DEFAULT,
88 buffer));
89
90 // Perform operation
91 HIPSPARSE_CHECK(hipsparseSpMV(handle,
92 transA,
93 &alpha,
94 matA,
95 vecX,
96 &beta,
97 vecY,
98 HIP_R_32F,
99 HIPSPARSE_MV_ALG_DEFAULT,
100 buffer));
101
102 // Copy device to host
103 HIP_CHECK(hipMemcpy(hy, dy, sizeof(float) * m, hipMemcpyDeviceToHost));
104
105 // Destroy matrix descriptors and handles
106 HIPSPARSE_CHECK(hipsparseDestroySpMat(matA));
107 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecX));
108 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
109 HIPSPARSE_CHECK(hipsparseDestroy(handle));
110
111 HIP_CHECK(hipFree(buffer));
112 HIP_CHECK(hipFree(dcsrRowPtr));
113 HIP_CHECK(hipFree(dcsrColInd));
114 HIP_CHECK(hipFree(dcsrVal));
115 HIP_CHECK(hipFree(dx));
116 HIP_CHECK(hipFree(dy));
117
118 return 0;
119}
hipsparseSpMM_bufferSize()#
-
hipsparseStatus_t hipsparseSpMM_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstDnMatDescr_t matB, const void *beta, const hipsparseDnMatDescr_t matC, hipDataType computeType, hipsparseSpMMAlg_t alg, size_t *pBufferSizeInBytes)#
hipsparseSpMM_bufferSizecomputes the required user allocated buffer size needed when computing the sparse matrix multiplication with a dense matrix:\[ C := \alpha \cdot op(A) \cdot op(B) + \beta \cdot C, \]where \(op(A)\) is a sparse \(m \times k\) matrix in CSR, COO, BSR or Blocked ELL storage format, \(B\) is a dense matrix of size \(k \times n\) and \(C\) is a dense matrix of size \(m \times n\).hipsparseSpMM_bufferSizesupports multiple combinations of data types and compute types. See hipsparseSpMM for a complete listing of all the data type and compute type combinations available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] matrix operation type.
opB – [in] matrix operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] sparse matrix descriptor.
matB – [in] dense matrix descriptor.
beta – [in] scalar \(\beta\).
matC – [in] dense matrix descriptor.
computeType – [in] floating point precision for the SpMM computation.
alg – [in] SpMM algorithm for the SpMM computation.
pBufferSizeInBytes – [out] number of bytes of the temporary storage buffer.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_NOT_INITIALIZED –
handleis not initialized.HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,matA,matB,matC,beta, orpBufferSizeInBytesis nullptr, oropAoropBis invalid.HIPSPARSE_STATUS_NOT_SUPPORTED –
opA,opB,computeTypeoralgis currently not supported.
hipsparseSpMM_preprocess()#
-
hipsparseStatus_t hipsparseSpMM_preprocess(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstDnMatDescr_t matB, const void *beta, const hipsparseDnMatDescr_t matC, hipDataType computeType, hipsparseSpMMAlg_t alg, void *externalBuffer)#
hipsparseSpMM_preprocessperforms the required preprocessing used when computing the sparse matrix multiplication with a dense matrix:\[ C := \alpha \cdot op(A) \cdot op(B) + \beta \cdot C, \]where \(op(A)\) is a sparse \(m \times k\) matrix in CSR, COO, BSR or Blocked ELL storage format, \(B\) is a dense matrix of size \(k \times n\) and \(C\) is a dense matrix of size \(m \times n\).hipsparseSpMM_preprocesssupports multiple combinations of data types and compute types. See hipsparseSpMM for a complete listing of all the data type and compute type combinations available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] matrix operation type.
opB – [in] matrix operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] matrix descriptor.
matB – [in] matrix descriptor.
beta – [in] scalar \(\beta\).
matC – [in] matrix descriptor.
computeType – [in] floating point precision for the SpMM computation.
alg – [in] SpMM algorithm for the SpMM computation.
externalBuffer – [out] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,matA,matB,matC,beta, orexternalBufferpointer is invalid.HIPSPARSE_STATUS_NOT_SUPPORTED –
opA,opB,computeTypeoralgis currently not supported.
hipsparseSpMM()#
-
hipsparseStatus_t hipsparseSpMM(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstDnMatDescr_t matB, const void *beta, const hipsparseDnMatDescr_t matC, hipDataType computeType, hipsparseSpMMAlg_t alg, void *externalBuffer)#
Compute the sparse matrix multiplication with a dense matrix.
hipsparseSpMMmultiplies the scalar \(\alpha\) with a sparse \(m \times k\) matrix \(op(A)\), defined in CSR, COO, BSR or Blocked ELL storage format, and the dense \(k \times n\) matrix \(op(B)\) and adds the result to the dense \(m \times n\) matrix \(C\) that is multiplied by the scalar \(\beta\), such that\[ C := \alpha \cdot op(A) \cdot op(B) + \beta \cdot C, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans_A == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ A^T, & \text{if trans_A == HIPSPARSE_OPERATION_TRANSPOSE} \\ A^H, & \text{if trans_A == HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE} \end{array} \right. \end{split}\]and\[\begin{split} op(B) = \left\{ \begin{array}{ll} B, & \text{if trans_B == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ B^T, & \text{if trans_B == HIPSPARSE_OPERATION_TRANSPOSE} \\ B^H, & \text{if trans_B == HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE} \end{array} \right. \end{split}\]Both \(B\) and \(C\) can be in row or column order.hipsparseSpMMrequires three stages to complete. First, the user calls hipsparseSpMM_bufferSize to determine the size of the required temporary storage buffer. Next, the user allocates this buffer and calls hipsparseSpMM_preprocess which will perform analysis on the sparse matrix \(op(A)\). Finally, the user callshipsparseSpMMto perform the actual computation. The buffer size and preprecess routines only need to be called once for a given sparse matrix \(op(A)\) while the computation routine can be repeatedly used with different \(B\) and \(C\) matrices. Once all calls tohipsparseSpMMare complete, the temporary buffer can be deallocated.As noted above, both \(B\) and \(C\) can be in row or column order (this includes mixing the order so that \(B\) is row order and \(C\) is column order and vice versa). For best performance, use row order for both \(B\) and \(C\) as this provides the best memory access.
hipsparseSpMMsupports multiple different algorithms. These algorithms have different trade offs depending on the sparsity pattern of the matrix, whether or not the results need to be deterministic, and how many times the sparse-matrix product will be performed.CSR Algorithms
HIPSPARSE_SPMM_CSR_ALG1
HIPSPARSE_SPMM_CSR_ALG2
HIPSPARSE_SPMM_CSR_ALG3
COO Algorithms
HIPSPARSE_SPMM_COO_ALG1
HIPSPARSE_SPMM_COO_ALG2
HIPSPARSE_SPMM_COO_ALG3
HIPSPARSE_SPMM_COO_ALG4
ELL Algorithms
HIPSPARSE_SPMM_BLOCKED_ELL_ALG1
BSR Algorithms
CUSPARSE_SPMM_BSR_ALG1
One can also pass HIPSPARSE_SPMM_ALG_DEFAULT which will automatically select from the algorithms listed above based on the sparse matrix format.
When A is transposed,
hipsparseSpMMwill revert to using HIPSPARSE_SPMM_CSR_ALG2 for CSR format and HIPSPARSE_SPMM_COO_ALG1 for COO format regardless of algorithm selected.hipsparseSpMMsupports multiple combinations of data types and compute types. The tables below indicate the currently supported different data types that can be used for for the sparse matrix \(op(A)\) and the dense matrices \(op(B)\) and \(C\) and the compute type for \(\alpha\) and \(\beta\). The advantage of using different data types is to save on memory bandwidth and storage when a user application allows while performing the actual computation in a higher precision.hipsparseSpMMsupports HIPSPARSE_INDEX_32I and HIPSPARSE_INDEX_64I index precisions for storing the row pointer and column indices arrays of the sparse matrices.- Uniform Precisions:
A / B / C / compute_type
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
- Mixed precisions:
A / B
C
compute_type
HIP_R_8I
HIP_R_32I
HIP_R_32I
HIP_R_8I
HIP_R_32F
HIP_R_32F
HIP_R_16F
HIP_R_32F
HIP_R_32F
HIP_R_16BF
HIP_R_32F
HIP_R_32F
hipsparseSpMMalso supports batched computation for CSR and COO matrices. There are three supported batch modes:\[\begin{split} C_i = A \times B_i \\ C_i = A_i \times B \\ C_i = A_i \times B_i \end{split}\]The batch mode is determined by the batch count and stride passed for each matrix. For example to use the first batch mode ( \(C_i = A \times B_i\)) with 100 batches for non-transposed \(A\), \(B\), and \(C\), one passes:
\[\begin{split} batchCountA=1 \\ batchCountB=100 \\ batchCountC=100 \\ offsetsBatchStrideA=0 \\ columnsValuesBatchStrideA=0 \\ batchStrideB=k*n \\ batchStrideC=m*n \end{split}\]To use the second batch mode ( \(C_i = A_i \times B\)) one could use:\[\begin{split} batchCountA=100 \\ batchCountB=1 \\ batchCountC=100 \\ offsetsBatchStrideA=m+1 \\ columnsValuesBatchStrideA=nnz \\ batchStrideB=0 \\ batchStrideC=m*n \end{split}\]And to use the third batch mode ( \(C_i = A_i \times B_i\)) one could use:\[\begin{split} batchCountA=100 \\ batchCountB=100 \\ batchCountC=100 \\ offsetsBatchStrideA=m+1 \\ columnsValuesBatchStrideA=nnz \\ batchStrideB=k*n \\ batchStrideC=m*n \end{split}\]See examples below.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] matrix operation type.
opB – [in] matrix operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] matrix descriptor.
matB – [in] matrix descriptor.
beta – [in] scalar \(\beta\).
matC – [in] matrix descriptor.
computeType – [in] floating point precision for the SpMM computation.
alg – [in] SpMM algorithm for the SpMM computation.
externalBuffer – [out] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,matA,matB,matC,beta, orexternalBufferpointer is invalid.HIPSPARSE_STATUS_NOT_SUPPORTED –
opA,opB,computeTypeoralgis currently not supported.
1int main(int argc, char* argv[])
2{
3 // A, B, and C are m×k, k×n, and m×n
4 int m = 3, n = 5, k = 4;
5 int ldb = n, ldc = n;
6 int nnz_A = 8, nnz_B = 20, nnz_C = 15;
7 hipsparseOperation_t transA = HIPSPARSE_OPERATION_NON_TRANSPOSE;
8 hipsparseOperation_t transB = HIPSPARSE_OPERATION_NON_TRANSPOSE;
9 hipsparseOperation_t transC = HIPSPARSE_OPERATION_NON_TRANSPOSE;
10 hipsparseOrder_t order = HIPSPARSE_ORDER_ROW;
11
12 // alpha and beta
13 float alpha = 0.5f;
14 float beta = 0.25f;
15
16 std::vector<int> hcsrRowPtr = {0, 3, 5, 8};
17 std::vector<int> hcsrColInd = {0, 1, 3, 1, 2, 0, 2, 3};
18 std::vector<float> hcsrVal = {1, 2, 3, 4, 5, 6, 7, 8};
19
20 std::vector<float> hB(nnz_B, 1.0f);
21 std::vector<float> hC(nnz_C, 1.0f);
22
23 int* dcsrRowPtr;
24 int* dcsrColInd;
25 float* dcsrVal;
26 HIP_CHECK(hipMalloc((void**)&dcsrRowPtr, sizeof(int) * (m + 1)));
27 HIP_CHECK(hipMalloc((void**)&dcsrColInd, sizeof(int) * nnz_A));
28 HIP_CHECK(hipMalloc((void**)&dcsrVal, sizeof(float) * nnz_A));
29
30 HIP_CHECK(
31 hipMemcpy(dcsrRowPtr, hcsrRowPtr.data(), sizeof(int) * (m + 1), hipMemcpyHostToDevice));
32 HIP_CHECK(hipMemcpy(dcsrColInd, hcsrColInd.data(), sizeof(int) * nnz_A, hipMemcpyHostToDevice));
33 HIP_CHECK(hipMemcpy(dcsrVal, hcsrVal.data(), sizeof(float) * nnz_A, hipMemcpyHostToDevice));
34
35 hipsparseHandle_t handle;
36 HIPSPARSE_CHECK(hipsparseCreate(&handle));
37
38 hipsparseSpMatDescr_t matA;
39 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
40 m,
41 k,
42 nnz_A,
43 dcsrRowPtr,
44 dcsrColInd,
45 dcsrVal,
46 HIPSPARSE_INDEX_32I,
47 HIPSPARSE_INDEX_32I,
48 HIPSPARSE_INDEX_BASE_ZERO,
49 HIP_R_32F));
50
51 // Allocate memory for the matrix B
52 float* dB;
53 HIP_CHECK(hipMalloc((void**)&dB, sizeof(float) * nnz_B));
54 HIP_CHECK(hipMemcpy(dB, hB.data(), sizeof(float) * nnz_B, hipMemcpyHostToDevice));
55
56 hipsparseDnMatDescr_t matB;
57 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matB, k, n, ldb, dB, HIP_R_32F, order));
58
59 // Allocate memory for the resulting matrix C
60 float* dC;
61 HIP_CHECK(hipMalloc((void**)&dC, sizeof(float) * nnz_C));
62 HIP_CHECK(hipMemcpy(dC, hC.data(), sizeof(float) * nnz_C, hipMemcpyHostToDevice));
63
64 hipsparseDnMatDescr_t matC;
65 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matC, m, n, ldc, dC, HIP_R_32F, HIPSPARSE_ORDER_ROW));
66
67 // Compute buffersize
68 size_t bufferSize;
69 HIPSPARSE_CHECK(hipsparseSpMM_bufferSize(handle,
70 transA,
71 transB,
72 &alpha,
73 matA,
74 matB,
75 &beta,
76 matC,
77 HIP_R_32F,
78 HIPSPARSE_MM_ALG_DEFAULT,
79 &bufferSize));
80
81 void* buffer;
82 HIP_CHECK(hipMalloc(&buffer, bufferSize));
83
84 // Preprocess operation (Optional)
85 HIPSPARSE_CHECK(hipsparseSpMM_preprocess(handle,
86 transA,
87 transB,
88 &alpha,
89 matA,
90 matB,
91 &beta,
92 matC,
93 HIP_R_32F,
94 HIPSPARSE_MM_ALG_DEFAULT,
95 buffer));
96
97 // Perform operation
98 HIPSPARSE_CHECK(hipsparseSpMM(handle,
99 transA,
100 transB,
101 &alpha,
102 matA,
103 matB,
104 &beta,
105 matC,
106 HIP_R_32F,
107 HIPSPARSE_MM_ALG_DEFAULT,
108 buffer));
109
110 // Copy device to host
111 HIP_CHECK(hipMemcpy(hC.data(), dC, sizeof(float) * nnz_C, hipMemcpyDeviceToHost));
112
113 // Destroy matrix descriptors and handles
114 HIPSPARSE_CHECK(hipsparseDestroySpMat(matA));
115 HIPSPARSE_CHECK(hipsparseDestroyDnMat(matB));
116 HIPSPARSE_CHECK(hipsparseDestroyDnMat(matC));
117 HIPSPARSE_CHECK(hipsparseDestroy(handle));
118
119 HIP_CHECK(hipFree(buffer));
120 HIP_CHECK(hipFree(dcsrRowPtr));
121 HIP_CHECK(hipFree(dcsrColInd));
122 HIP_CHECK(hipFree(dcsrVal));
123 HIP_CHECK(hipFree(dB));
124 HIP_CHECK(hipFree(dC));
125
126 return 0;
127}
1int main(int argc, char* argv[])
2{
3 // A, B, and C are m×k, k×n, and m×n
4 int m = 3, n = 5, k = 4;
5 int ldb = n, ldc = n;
6 int nnz_A = 8, nnz_B = 20, nnz_C = 15;
7 hipsparseOperation_t transA = HIPSPARSE_OPERATION_NON_TRANSPOSE;
8 hipsparseOperation_t transB = HIPSPARSE_OPERATION_NON_TRANSPOSE;
9 hipsparseOperation_t transC = HIPSPARSE_OPERATION_NON_TRANSPOSE;
10 hipsparseOrder_t order = HIPSPARSE_ORDER_ROW;
11
12 // alpha and beta
13 float alpha = 0.5;
14 float beta = 0.25;
15
16 int hcsrRowPtr[] = {0, 3, 5, 8};
17 int hcsrColInd[] = {0, 1, 3, 1, 2, 0, 2, 3};
18 float hcsrVal[] = {1, 2, 3, 4, 5, 6, 7, 8};
19
20 float* hB = (float*)malloc(nnz_B * sizeof(float));
21 float* hC = (float*)malloc(nnz_C * sizeof(float));
22
23 int* dcsrRowPtr;
24 int* dcsrColInd;
25 float* dcsrVal;
26 HIP_CHECK(hipMalloc((void**)&dcsrRowPtr, sizeof(int) * (m + 1)));
27 HIP_CHECK(hipMalloc((void**)&dcsrColInd, sizeof(int) * nnz_A));
28 HIP_CHECK(hipMalloc((void**)&dcsrVal, sizeof(float) * nnz_A));
29
30 HIP_CHECK(hipMemcpy(dcsrRowPtr, hcsrRowPtr, sizeof(int) * (m + 1), hipMemcpyHostToDevice));
31 HIP_CHECK(hipMemcpy(dcsrColInd, hcsrColInd, sizeof(int) * nnz_A, hipMemcpyHostToDevice));
32 HIP_CHECK(hipMemcpy(dcsrVal, hcsrVal, sizeof(float) * nnz_A, hipMemcpyHostToDevice));
33
34 hipsparseHandle_t handle;
35 HIPSPARSE_CHECK(hipsparseCreate(&handle));
36
37 hipsparseSpMatDescr_t matA;
38 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
39 m,
40 k,
41 nnz_A,
42 dcsrRowPtr,
43 dcsrColInd,
44 dcsrVal,
45 HIPSPARSE_INDEX_32I,
46 HIPSPARSE_INDEX_32I,
47 HIPSPARSE_INDEX_BASE_ZERO,
48 HIP_R_32F));
49
50 // Allocate memory for the matrix B
51 float* dB;
52 HIP_CHECK(hipMalloc((void**)&dB, sizeof(float) * nnz_B));
53 HIP_CHECK(hipMemcpy(dB, hB, sizeof(float) * nnz_B, hipMemcpyHostToDevice));
54
55 hipsparseDnMatDescr_t matB;
56 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matB, k, n, ldb, dB, HIP_R_32F, order));
57
58 // Allocate memory for the resulting matrix C
59 float* dC;
60 HIP_CHECK(hipMalloc((void**)&dC, sizeof(float) * nnz_C));
61 HIP_CHECK(hipMemcpy(dC, hC, sizeof(float) * nnz_C, hipMemcpyHostToDevice));
62
63 hipsparseDnMatDescr_t matC;
64 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matC, m, n, ldc, dC, HIP_R_32F, HIPSPARSE_ORDER_ROW));
65
66 // Compute buffersize
67 size_t bufferSize;
68 HIPSPARSE_CHECK(hipsparseSpMM_bufferSize(handle,
69 transA,
70 transB,
71 &alpha,
72 matA,
73 matB,
74 &beta,
75 matC,
76 HIP_R_32F,
77 HIPSPARSE_MM_ALG_DEFAULT,
78 &bufferSize));
79
80 void* buffer;
81 HIP_CHECK(hipMalloc(&buffer, bufferSize));
82
83 // Preprocess operation (Optional)
84 HIPSPARSE_CHECK(hipsparseSpMM_preprocess(handle,
85 transA,
86 transB,
87 &alpha,
88 matA,
89 matB,
90 &beta,
91 matC,
92 HIP_R_32F,
93 HIPSPARSE_MM_ALG_DEFAULT,
94 buffer));
95
96 // Perform operation
97 HIPSPARSE_CHECK(hipsparseSpMM(handle,
98 transA,
99 transB,
100 &alpha,
101 matA,
102 matB,
103 &beta,
104 matC,
105 HIP_R_32F,
106 HIPSPARSE_MM_ALG_DEFAULT,
107 buffer));
108
109 // Copy device to host
110 HIP_CHECK(hipMemcpy(hC, dC, sizeof(float) * nnz_C, hipMemcpyDeviceToHost));
111
112 // Destroy matrix descriptors and handles
113 HIPSPARSE_CHECK(hipsparseDestroySpMat(matA));
114 HIPSPARSE_CHECK(hipsparseDestroyDnMat(matB));
115 HIPSPARSE_CHECK(hipsparseDestroyDnMat(matC));
116 HIPSPARSE_CHECK(hipsparseDestroy(handle));
117
118 HIP_CHECK(hipFree(buffer));
119 HIP_CHECK(hipFree(dcsrRowPtr));
120 HIP_CHECK(hipFree(dcsrColInd));
121 HIP_CHECK(hipFree(dcsrVal));
122 HIP_CHECK(hipFree(dB));
123 HIP_CHECK(hipFree(dC));
124
125 return 0;
126}
hipsparseSpGEMM_createDescr()#
-
hipsparseStatus_t hipsparseSpGEMM_createDescr(hipsparseSpGEMMDescr_t *descr)#
hipsparseSpGEMM_createDescrcreates a sparse matrix sparse matrix product descriptor. It should be destroyed at the end using hipsparseSpGEMM_destroyDescr().
hipsparseSpGEMM_destroyDescr()#
-
hipsparseStatus_t hipsparseSpGEMM_destroyDescr(hipsparseSpGEMMDescr_t descr)#
hipsparseSpGEMM_destroyDescrdestroys a sparse matrix sparse matrix product descriptor and releases all resources used by the descriptor.
hipsparseSpGEMM_workEstimation()#
-
hipsparseStatus_t hipsparseSpGEMM_workEstimation(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstSpMatDescr_t matB, const void *beta, hipsparseSpMatDescr_t matC, hipDataType computeType, hipsparseSpGEMMAlg_t alg, hipsparseSpGEMMDescr_t spgemmDescr, size_t *bufferSize1, void *externalBuffer1)#
Work estimation step of the sparse matrix sparse matrix product:
\[ C' := \alpha \cdot op(A) \cdot op(B) + \beta \cdot C, \]where \(C'\), \(A\), \(B\), \(C\) are sparse matrices and \(C'\) and \(C\) have the same sparsity pattern.hipsparseSpGEMM_workEstimationis called twice. We call it to compute the size of the first required user allocated buffer. After this buffer size is determined, the user allocates it and callshipsparseSpGEMM_workEstimationa second time with the newly allocated buffer passed in. This second call inspects the matrices \(A\) and \(B\) to determine the number of intermediate products that will result from multipltying \(A\) and \(B\) together.hipsparseSpGEMM_workEstimationsupports multiple combinations of data types and compute types. See hipsparseSpGEMM_copy for a complete listing of all the data type and compute type combinations available.- Example (See full example below)
void* dBuffer1 = NULL; size_t bufferSize1 = 0; hipsparseSpGEMMDescr_t spgemmDesc; hipsparseSpGEMM_createDescr(&spgemmDesc); size_t bufferSize1 = 0; hipsparseSpGEMM_workEstimation(handle, opA, opB, &alpha, matA, matB, &beta, matC, computeType, HIPSPARSE_SPGEMM_DEFAULT, spgemmDesc, &bufferSize1, NULL); hipMalloc((void**) &dBuffer1, bufferSize1); // Determine number of intermediate product when computing A * B hipsparseSpGEMM_workEstimation(handle, opA, opB, &alpha, matA, matB, &beta, matC, computeType, HIPSPARSE_SPGEMM_DEFAULT, spgemmDesc, &bufferSize1, dBuffer1);
- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] sparse matrix \(A\) operation type.
opB – [in] sparse matrix \(B\) operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] sparse matrix \(A\) descriptor.
matB – [in] sparse matrix \(B\) descriptor.
beta – [in] scalar \(\beta\).
matC – [out] sparse matrix \(C\) descriptor.
computeType – [in] floating point precision for the SpGEMM computation.
alg – [in] SpGEMM algorithm for the SpGEMM computation.
spgemmDescr – [in] SpGEMM descriptor.
bufferSize1 – [out] number of bytes of the temporary storage buffer.
externalBuffer1 – [in] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_NOT_INITIALIZED –
handleis not initialized.HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,beta,matA,matB,matC,spgemmDescrorbufferSize1is nullptr, oropAoropBis invalid.HIPSPARSE_STATUS_ALLOC_FAILED – additional buffer for long rows could not be allocated.
HIPSPARSE_STATUS_NOT_SUPPORTED –
opA!= HIPSPARSE_OPERATION_NON_TRANSPOSE oropB!= HIPSPARSE_OPERATION_NON_TRANSPOSE.
hipsparseSpGEMM_compute()#
-
hipsparseStatus_t hipsparseSpGEMM_compute(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstSpMatDescr_t matB, const void *beta, hipsparseSpMatDescr_t matC, hipDataType computeType, hipsparseSpGEMMAlg_t alg, hipsparseSpGEMMDescr_t spgemmDescr, size_t *bufferSize2, void *externalBuffer2)#
Compute step of the sparse matrix sparse matrix product:
\[ C' := \alpha \cdot op(A) \cdot op(B) + \beta \cdot C, \]where \(C'\), \(A\), \(B\), \(C\) are sparse matrices and \(C'\) and \(C\) have the same sparsity pattern.hipsparseSpGEMM_computeis called twice. First to compute the size of the second required user allocated buffer. After this buffer size is determined, the user allocates it and callshipsparseSpGEMM_computea second time with the newly allocated buffer passed in. This second call performs the actual computation of \(C' = \alpha \cdot A \cdot B\) (the result is stored in the temporary buffers).hipsparseSpGEMM_computesupports multiple combinations of data types and compute types. See hipsparseSpGEMM_copy for a complete listing of all the data type and compute type combinations available.- Example (See full example below)
void* dBuffer2 = NULL; size_t bufferSize2 = 0; size_t bufferSize2 = 0; hipsparseSpGEMM_compute(handle, opA, opB, &alpha, matA, matB, &beta, matC, computeType, HIPSPARSE_SPGEMM_DEFAULT, spgemmDesc, &bufferSize2, NULL); hipMalloc((void**) &dBuffer2, bufferSize2); // compute the intermediate product of A * B hipsparseSpGEMM_compute(handle, opA, opB, &alpha, matA, matB, &beta, matC, computeType, HIPSPARSE_SPGEMM_DEFAULT, spgemmDesc, &bufferSize2, dBuffer2);
- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] sparse matrix \(A\) operation type.
opB – [in] sparse matrix \(B\) operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] sparse matrix \(A\) descriptor.
matB – [in] sparse matrix \(B\) descriptor.
beta – [in] scalar \(\beta\).
matC – [out] sparse matrix \(C\) descriptor.
computeType – [in] floating point precision for the SpGEMM computation.
alg – [in] SpGEMM algorithm for the SpGEMM computation.
spgemmDescr – [in] SpGEMM descriptor.
bufferSize2 – [out] number of bytes of the temporary storage buffer.
externalBuffer2 – [in] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,beta,matA,matB,matCorbufferSize2pointer is invalid.HIPSPARSE_STATUS_ALLOC_FAILED – additional buffer for long rows could not be allocated.
HIPSPARSE_STATUS_NOT_SUPPORTED –
opA!= HIPSPARSE_OPERATION_NON_TRANSPOSE oropB!= HIPSPARSE_OPERATION_NON_TRANSPOSE.
hipsparseSpGEMM_copy()#
-
hipsparseStatus_t hipsparseSpGEMM_copy(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstSpMatDescr_t matB, const void *beta, hipsparseSpMatDescr_t matC, hipDataType computeType, hipsparseSpGEMMAlg_t alg, hipsparseSpGEMMDescr_t spgemmDescr)#
Copy step of the sparse matrix sparse matrix product:
\[ C' := \alpha \cdot op(A) \cdot op(B) + \beta \cdot C, \]where \(C'\), \(A\), \(B\), \(C\) are sparse matrices and \(C'\) and \(C\) have the same sparsity pattern.hipsparseSpGEMM_copyis called once to copy the results (that are currently stored in the temporary arrays) to the output sparse matrix. If \(\beta != 0\), then the \(beta \cdot C\) portion of the computation: \(C' = \alpha \cdot A \cdot B + \beta * C\) is handled. This is possible because \(C'\) and \(C\) must have the same sparsity pattern.hipsparseSpGEMM_copysupports multiple combinations of data types and compute types. The tables below indicate the currently supported different data types that can be used for for the sparse matrices \(op(A)\), \(op(B)\), \(C\), and \(C'\) and the compute type for \(\alpha\) and \(\beta\). The advantage of using different data types is to save on memory bandwidth and storage when a user application allows while performing the actual computation in a higher precision.hipsparseSpGEMM_copysupports HIPSPARSE_INDEX_32I and HIPSPARSE_INDEX_64I index precisions for storing the row pointer and row/column indices arrays of the sparse matrices.- Uniform Precisions:
A / B / C / C’ / compute_type
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
Note
The two user allocated temporary buffers can only be freed after the call to
hipsparseSpGEMM_copy- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] sparse matrix \(A\) operation type.
opB – [in] sparse matrix \(B\) operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] sparse matrix \(A\) descriptor.
matB – [in] sparse matrix \(B\) descriptor.
beta – [in] scalar \(\beta\).
matC – [out] sparse matrix \(C\) descriptor.
computeType – [in] floating point precision for the SpGEMM computation.
alg – [in] SpGEMM algorithm for the SpGEMM computation.
spgemmDescr – [in] SpGEMM descriptor.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,beta,matA,matB,matCpointer is invalid.HIPSPARSE_STATUS_ALLOC_FAILED – additional buffer for long rows could not be allocated.
HIPSPARSE_STATUS_NOT_SUPPORTED –
opA!= HIPSPARSE_OPERATION_NON_TRANSPOSE oropB!= HIPSPARSE_OPERATION_NON_TRANSPOSE.
1int main(int argc, char* argv[])
2{
3 int m = 2;
4 int k = 2;
5 int n = 3;
6 int nnzA = 4;
7 int nnzB = 4;
8
9 float alpha{1.0f};
10 float beta{0.0f};
11
12 hipsparseOperation_t opA = HIPSPARSE_OPERATION_NON_TRANSPOSE;
13 hipsparseOperation_t opB = HIPSPARSE_OPERATION_NON_TRANSPOSE;
14 hipDataType computeType = HIP_R_32F;
15
16 // A, B, and C are m×k, k×n, and m×n
17
18 // A
19 std::vector<int> hcsrRowPtrA = {0, 2, 4};
20 std::vector<int> hcsrColIndA = {0, 1, 0, 1};
21 std::vector<float> hcsrValA = {1.0f, 2.0f, 3.0f, 4.0f};
22
23 // B
24 std::vector<int> hcsrRowPtrB = {0, 2, 4};
25 std::vector<int> hcsrColIndB = {1, 2, 0, 2};
26 std::vector<float> hcsrValB = {5.0f, 6.0f, 7.0f, 8.0f};
27
28 // Device memory management: Allocate and copy A, B
29 int* dcsrRowPtrA;
30 int* dcsrColIndA;
31 float* dcsrValA;
32 int* dcsrRowPtrB;
33 int* dcsrColIndB;
34 float* dcsrValB;
35 int* dcsrRowPtrC;
36 int* dcsrColIndC;
37 float* dcsrValC;
38 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrA, (m + 1) * sizeof(int)));
39 HIP_CHECK(hipMalloc((void**)&dcsrColIndA, nnzA * sizeof(int)));
40 HIP_CHECK(hipMalloc((void**)&dcsrValA, nnzA * sizeof(float)));
41 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrB, (k + 1) * sizeof(int)));
42 HIP_CHECK(hipMalloc((void**)&dcsrColIndB, nnzB * sizeof(int)));
43 HIP_CHECK(hipMalloc((void**)&dcsrValB, nnzB * sizeof(float)));
44 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrC, (m + 1) * sizeof(int)));
45
46 HIP_CHECK(
47 hipMemcpy(dcsrRowPtrA, hcsrRowPtrA.data(), (m + 1) * sizeof(int), hipMemcpyHostToDevice));
48 HIP_CHECK(
49 hipMemcpy(dcsrColIndA, hcsrColIndA.data(), nnzA * sizeof(int), hipMemcpyHostToDevice));
50 HIP_CHECK(hipMemcpy(dcsrValA, hcsrValA.data(), nnzA * sizeof(float), hipMemcpyHostToDevice));
51 HIP_CHECK(
52 hipMemcpy(dcsrRowPtrB, hcsrRowPtrB.data(), (k + 1) * sizeof(int), hipMemcpyHostToDevice));
53 HIP_CHECK(
54 hipMemcpy(dcsrColIndB, hcsrColIndB.data(), nnzB * sizeof(int), hipMemcpyHostToDevice));
55 HIP_CHECK(hipMemcpy(dcsrValB, hcsrValB.data(), nnzB * sizeof(float), hipMemcpyHostToDevice));
56
57 hipsparseSpMatDescr_t matA, matB, matC;
58
59 hipsparseHandle_t handle;
60 HIPSPARSE_CHECK(hipsparseCreate(&handle));
61
62 // Create sparse matrix A in CSR format
63 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
64 m,
65 k,
66 nnzA,
67 dcsrRowPtrA,
68 dcsrColIndA,
69 dcsrValA,
70 HIPSPARSE_INDEX_32I,
71 HIPSPARSE_INDEX_32I,
72 HIPSPARSE_INDEX_BASE_ZERO,
73 HIP_R_32F));
74 HIPSPARSE_CHECK(hipsparseCreateCsr(&matB,
75 k,
76 n,
77 nnzB,
78 dcsrRowPtrB,
79 dcsrColIndB,
80 dcsrValB,
81 HIPSPARSE_INDEX_32I,
82 HIPSPARSE_INDEX_32I,
83 HIPSPARSE_INDEX_BASE_ZERO,
84 HIP_R_32F));
85 HIPSPARSE_CHECK(hipsparseCreateCsr(&matC,
86 m,
87 n,
88 0,
89 dcsrRowPtrC,
90 NULL,
91 NULL,
92 HIPSPARSE_INDEX_32I,
93 HIPSPARSE_INDEX_32I,
94 HIPSPARSE_INDEX_BASE_ZERO,
95 HIP_R_32F));
96
97 hipsparseSpGEMMDescr_t spgemmDesc;
98 HIPSPARSE_CHECK(hipsparseSpGEMM_createDescr(&spgemmDesc));
99
100 void* dBuffer1 = NULL;
101 void* dBuffer2 = NULL;
102 size_t bufferSize1 = 0;
103 size_t bufferSize2 = 0;
104
105 // Determine size of first user allocated buffer
106 HIPSPARSE_CHECK(hipsparseSpGEMM_workEstimation(handle,
107 opA,
108 opB,
109 &alpha,
110 matA,
111 matB,
112 &beta,
113 matC,
114 computeType,
115 HIPSPARSE_SPGEMM_DEFAULT,
116 spgemmDesc,
117 &bufferSize1,
118 NULL));
119 HIP_CHECK(hipMalloc((void**)&dBuffer1, bufferSize1));
120
121 // Inspect the matrices A and B to determine the number of intermediate product in
122 // C = alpha * A * B
123 HIPSPARSE_CHECK(hipsparseSpGEMM_workEstimation(handle,
124 opA,
125 opB,
126 &alpha,
127 matA,
128 matB,
129 &beta,
130 matC,
131 computeType,
132 HIPSPARSE_SPGEMM_DEFAULT,
133 spgemmDesc,
134 &bufferSize1,
135 dBuffer1));
136
137 // Determine size of second user allocated buffer
138 HIPSPARSE_CHECK(hipsparseSpGEMM_compute(handle,
139 opA,
140 opB,
141 &alpha,
142 matA,
143 matB,
144 &beta,
145 matC,
146 computeType,
147 HIPSPARSE_SPGEMM_DEFAULT,
148 spgemmDesc,
149 &bufferSize2,
150 NULL));
151 HIP_CHECK(hipMalloc((void**)&dBuffer2, bufferSize2));
152
153 // Compute C = alpha * A * B and store result in temporary buffers
154 HIPSPARSE_CHECK(hipsparseSpGEMM_compute(handle,
155 opA,
156 opB,
157 &alpha,
158 matA,
159 matB,
160 &beta,
161 matC,
162 computeType,
163 HIPSPARSE_SPGEMM_DEFAULT,
164 spgemmDesc,
165 &bufferSize2,
166 dBuffer2));
167
168 // Get matrix C non-zero entries C_nnz1
169 int64_t C_num_rows1, C_num_cols1, C_nnz1;
170 HIPSPARSE_CHECK(hipsparseSpMatGetSize(matC, &C_num_rows1, &C_num_cols1, &C_nnz1));
171
172 // Allocate the CSR structures for the matrix C
173 HIP_CHECK(hipMalloc((void**)&dcsrColIndC, C_nnz1 * sizeof(int)));
174 HIP_CHECK(hipMalloc((void**)&dcsrValC, C_nnz1 * sizeof(float)));
175
176 // Update matC with the new pointers
177 HIPSPARSE_CHECK(hipsparseCsrSetPointers(matC, dcsrRowPtrC, dcsrColIndC, dcsrValC));
178
179 // Copy the final products to the matrix C
180 HIPSPARSE_CHECK(hipsparseSpGEMM_copy(handle,
181 opA,
182 opB,
183 &alpha,
184 matA,
185 matB,
186 &beta,
187 matC,
188 computeType,
189 HIPSPARSE_SPGEMM_DEFAULT,
190 spgemmDesc));
191
192 std::vector<int> hcsrRowPtrC(m + 1);
193 std::vector<int> hcsrColIndC(C_nnz1);
194 std::vector<float> hcsrValC(C_nnz1);
195
196 // Copy back to the host
197 HIP_CHECK(
198 hipMemcpy(hcsrRowPtrC.data(), dcsrRowPtrC, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
199 HIP_CHECK(
200 hipMemcpy(hcsrColIndC.data(), dcsrColIndC, sizeof(int) * C_nnz1, hipMemcpyDeviceToHost));
201 HIP_CHECK(hipMemcpy(hcsrValC.data(), dcsrValC, sizeof(float) * C_nnz1, hipMemcpyDeviceToHost));
202
203 std::cout << "C" << std::endl;
204 for(int i = 0; i < m; i++)
205 {
206 int start = hcsrRowPtrC[i];
207 int end = hcsrRowPtrC[i + 1];
208
209 std::vector<float> temp(n, 0.0f);
210 for(int j = start; j < end; j++)
211 {
212 temp[hcsrColIndC[j]] = hcsrValC[j];
213 }
214
215 for(int j = 0; j < n; j++)
216 {
217 std::cout << temp[j] << " ";
218 }
219 std::cout << std::endl;
220 }
221 std::cout << std::endl;
222
223 // Destroy matrix descriptors and handles
224 HIPSPARSE_CHECK(hipsparseSpGEMM_destroyDescr(spgemmDesc));
225 HIPSPARSE_CHECK(hipsparseDestroySpMat(matA));
226 HIPSPARSE_CHECK(hipsparseDestroySpMat(matB));
227 HIPSPARSE_CHECK(hipsparseDestroySpMat(matC));
228 HIPSPARSE_CHECK(hipsparseDestroy(handle));
229
230 // Free device memory
231 HIP_CHECK(hipFree(dBuffer1));
232 HIP_CHECK(hipFree(dBuffer2));
233
234 return 0;
235}
1int main(int argc, char* argv[])
2{
3 int m = 2;
4 int k = 2;
5 int n = 3;
6 int nnzA = 4;
7 int nnzB = 4;
8
9 float alpha = 1.0;
10 float beta = 0.0;
11
12 hipsparseOperation_t opA = HIPSPARSE_OPERATION_NON_TRANSPOSE;
13 hipsparseOperation_t opB = HIPSPARSE_OPERATION_NON_TRANSPOSE;
14 hipDataType computeType = HIP_R_32F;
15
16 // A, B, and C are m×k, k×n, and m×n
17
18 // A
19 int hcsrRowPtrA[] = {0, 2, 4};
20 int hcsrColIndA[] = {0, 1, 0, 1};
21 float hcsrValA[] = {1.0, 2.0, 3.0, 4.0};
22
23 // B
24 int hcsrRowPtrB[] = {0, 2, 4};
25 int hcsrColIndB[] = {1, 2, 0, 2};
26 float hcsrValB[] = {5.0, 6.0, 7.0, 8.0};
27
28 // Device memory management: Allocate and copy A, B
29 int* dcsrRowPtrA;
30 int* dcsrColIndA;
31 float* dcsrValA;
32 int* dcsrRowPtrB;
33 int* dcsrColIndB;
34 float* dcsrValB;
35 int* dcsrRowPtrC;
36 int* dcsrColIndC;
37 float* dcsrValC;
38 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrA, (m + 1) * sizeof(int)));
39 HIP_CHECK(hipMalloc((void**)&dcsrColIndA, nnzA * sizeof(int)));
40 HIP_CHECK(hipMalloc((void**)&dcsrValA, nnzA * sizeof(float)));
41 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrB, (k + 1) * sizeof(int)));
42 HIP_CHECK(hipMalloc((void**)&dcsrColIndB, nnzB * sizeof(int)));
43 HIP_CHECK(hipMalloc((void**)&dcsrValB, nnzB * sizeof(float)));
44 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrC, (m + 1) * sizeof(int)));
45
46 HIP_CHECK(hipMemcpy(dcsrRowPtrA, hcsrRowPtrA, (m + 1) * sizeof(int), hipMemcpyHostToDevice));
47 HIP_CHECK(hipMemcpy(dcsrColIndA, hcsrColIndA, nnzA * sizeof(int), hipMemcpyHostToDevice));
48 HIP_CHECK(hipMemcpy(dcsrValA, hcsrValA, nnzA * sizeof(float), hipMemcpyHostToDevice));
49 HIP_CHECK(hipMemcpy(dcsrRowPtrB, hcsrRowPtrB, (k + 1) * sizeof(int), hipMemcpyHostToDevice));
50 HIP_CHECK(hipMemcpy(dcsrColIndB, hcsrColIndB, nnzB * sizeof(int), hipMemcpyHostToDevice));
51 HIP_CHECK(hipMemcpy(dcsrValB, hcsrValB, nnzB * sizeof(float), hipMemcpyHostToDevice));
52
53 hipsparseSpMatDescr_t matA, matB, matC;
54
55 hipsparseHandle_t handle;
56 HIPSPARSE_CHECK(hipsparseCreate(&handle));
57
58 // Create sparse matrix A in CSR format
59 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
60 m,
61 k,
62 nnzA,
63 dcsrRowPtrA,
64 dcsrColIndA,
65 dcsrValA,
66 HIPSPARSE_INDEX_32I,
67 HIPSPARSE_INDEX_32I,
68 HIPSPARSE_INDEX_BASE_ZERO,
69 HIP_R_32F));
70 HIPSPARSE_CHECK(hipsparseCreateCsr(&matB,
71 k,
72 n,
73 nnzB,
74 dcsrRowPtrB,
75 dcsrColIndB,
76 dcsrValB,
77 HIPSPARSE_INDEX_32I,
78 HIPSPARSE_INDEX_32I,
79 HIPSPARSE_INDEX_BASE_ZERO,
80 HIP_R_32F));
81 HIPSPARSE_CHECK(hipsparseCreateCsr(&matC,
82 m,
83 n,
84 0,
85 dcsrRowPtrC,
86 NULL,
87 NULL,
88 HIPSPARSE_INDEX_32I,
89 HIPSPARSE_INDEX_32I,
90 HIPSPARSE_INDEX_BASE_ZERO,
91 HIP_R_32F));
92
93 hipsparseSpGEMMDescr_t spgemmDesc;
94 HIPSPARSE_CHECK(hipsparseSpGEMM_createDescr(&spgemmDesc));
95
96 void* dBuffer1 = NULL;
97 void* dBuffer2 = NULL;
98 size_t bufferSize1 = 0;
99 size_t bufferSize2 = 0;
100
101 // Determine size of first user allocated buffer
102 HIPSPARSE_CHECK(hipsparseSpGEMM_workEstimation(handle,
103 opA,
104 opB,
105 &alpha,
106 matA,
107 matB,
108 &beta,
109 matC,
110 computeType,
111 HIPSPARSE_SPGEMM_DEFAULT,
112 spgemmDesc,
113 &bufferSize1,
114 NULL));
115 HIP_CHECK(hipMalloc((void**)&dBuffer1, bufferSize1));
116
117 // Inspect the matrices A and B to determine the number of intermediate product in
118 // C = alpha * A * B
119 HIPSPARSE_CHECK(hipsparseSpGEMM_workEstimation(handle,
120 opA,
121 opB,
122 &alpha,
123 matA,
124 matB,
125 &beta,
126 matC,
127 computeType,
128 HIPSPARSE_SPGEMM_DEFAULT,
129 spgemmDesc,
130 &bufferSize1,
131 dBuffer1));
132
133 // Determine size of second user allocated buffer
134 HIPSPARSE_CHECK(hipsparseSpGEMM_compute(handle,
135 opA,
136 opB,
137 &alpha,
138 matA,
139 matB,
140 &beta,
141 matC,
142 computeType,
143 HIPSPARSE_SPGEMM_DEFAULT,
144 spgemmDesc,
145 &bufferSize2,
146 NULL));
147 HIP_CHECK(hipMalloc((void**)&dBuffer2, bufferSize2));
148
149 // Compute C = alpha * A * B and store result in temporary buffers
150 HIPSPARSE_CHECK(hipsparseSpGEMM_compute(handle,
151 opA,
152 opB,
153 &alpha,
154 matA,
155 matB,
156 &beta,
157 matC,
158 computeType,
159 HIPSPARSE_SPGEMM_DEFAULT,
160 spgemmDesc,
161 &bufferSize2,
162 dBuffer2));
163
164 // Get matrix C non-zero entries C_nnz1
165 int64_t C_num_rows1, C_num_cols1, C_nnz1;
166 HIPSPARSE_CHECK(hipsparseSpMatGetSize(matC, &C_num_rows1, &C_num_cols1, &C_nnz1));
167
168 // Allocate the CSR structures for the matrix C
169 HIP_CHECK(hipMalloc((void**)&dcsrColIndC, C_nnz1 * sizeof(int)));
170 HIP_CHECK(hipMalloc((void**)&dcsrValC, C_nnz1 * sizeof(float)));
171
172 // Update matC with the new pointers
173 HIPSPARSE_CHECK(hipsparseCsrSetPointers(matC, dcsrRowPtrC, dcsrColIndC, dcsrValC));
174
175 // Copy the final products to the matrix C
176 HIPSPARSE_CHECK(hipsparseSpGEMM_copy(handle,
177 opA,
178 opB,
179 &alpha,
180 matA,
181 matB,
182 &beta,
183 matC,
184 computeType,
185 HIPSPARSE_SPGEMM_DEFAULT,
186 spgemmDesc));
187
188 int* hcsrRowPtrC = (int*)malloc((m + 1) * sizeof(int));
189 int* hcsrColIndC = (int*)malloc((C_nnz1) * sizeof(int));
190 float hcsrValC[C_nnz1];
191
192 // Copy back to the host
193 HIP_CHECK(hipMemcpy(hcsrRowPtrC, dcsrRowPtrC, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
194 HIP_CHECK(hipMemcpy(hcsrColIndC, dcsrColIndC, sizeof(int) * C_nnz1, hipMemcpyDeviceToHost));
195 HIP_CHECK(hipMemcpy(hcsrValC, dcsrValC, sizeof(float) * C_nnz1, hipMemcpyDeviceToHost));
196
197 printf("C\n");
198 for(int i = 0; i < m; i++)
199 {
200 int start = hcsrRowPtrC[i];
201 int end = hcsrRowPtrC[i + 1];
202
203 float* temp = (float*)malloc(n * sizeof(float));
204 for(int j = start; j < end; j++)
205 {
206 temp[hcsrColIndC[j]] = hcsrValC[j];
207 }
208
209 for(int j = 0; j < n; j++)
210 {
211 printf("%f ", temp[j]);
212 }
213 printf("\n");
214 }
215 printf("\n");
216
217 // Destroy matrix descriptors and handles
218 HIPSPARSE_CHECK(hipsparseSpGEMM_destroyDescr(spgemmDesc));
219 HIPSPARSE_CHECK(hipsparseDestroySpMat(matA));
220 HIPSPARSE_CHECK(hipsparseDestroySpMat(matB));
221 HIPSPARSE_CHECK(hipsparseDestroySpMat(matC));
222 HIPSPARSE_CHECK(hipsparseDestroy(handle));
223
224 // Free device memory
225 HIP_CHECK(hipFree(dBuffer1));
226 HIP_CHECK(hipFree(dBuffer2));
227
228 return 0;
229}
hipsparseSpGEMMreuse_workEstimation()#
-
hipsparseStatus_t hipsparseSpGEMMreuse_workEstimation(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, hipsparseConstSpMatDescr_t matA, hipsparseConstSpMatDescr_t matB, hipsparseSpMatDescr_t matC, hipsparseSpGEMMAlg_t alg, hipsparseSpGEMMDescr_t spgemmDescr, size_t *bufferSize1, void *externalBuffer1)#
Work estimation step of the sparse matrix sparse matrix product:
\[ C' := \alpha \cdot op(A) \cdot op(B) + \beta \cdot C, \]where \(C'\), \(A\), \(B\), \(C\) are sparse matrices and \(C'\) and \(C\) have the same sparsity pattern.hipsparseSpGEMMreuse_workEstimationis called twice. We call it to compute the size of the first required user allocated buffer. After this buffer size is determined, the user allocates it and callshipsparseSpGEMMreuse_workEstimationa second time with the newly allocated buffer passed in. This second call inspects the matrices \(A\) and \(B\) to determine the number of intermediate products that will result from multipltying \(A\) and \(B\) together.- Example (See full example below)
void* dBuffer1 = NULL; size_t bufferSize1 = 0; hipsparseSpGEMMDescr_t spgemmDesc; hipsparseSpGEMM_createDescr(&spgemmDesc); size_t bufferSize1 = 0; hipsparseSpGEMMreuse_workEstimation(handle, opA, opB, matA, matB, matC, HIPSPARSE_SPGEMM_DEFAULT, spgemmDesc, &bufferSize1, NULL); hipMalloc((void**) &dBuffer1, bufferSize1); // Determine number of intermediate product when computing A * B hipsparseSpGEMMreuse_workEstimation(handle, opA, opB, matA, matB, matC, HIPSPARSE_SPGEMM_DEFAULT, spgemmDesc, &bufferSize1, dBuffer1);
- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] sparse matrix \(A\) operation type.
opB – [in] sparse matrix \(B\) operation type.
matA – [in] sparse matrix \(A\) descriptor.
matB – [in] sparse matrix \(B\) descriptor.
matC – [out] sparse matrix \(C\) descriptor.
alg – [in] SpGEMM algorithm for the SpGEMM computation.
spgemmDescr – [in] SpGEMM descriptor.
bufferSize1 – [out] number of bytes of the temporary storage buffer.
externalBuffer1 – [in] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,matA,matB,matCorbufferSize1pointer is invalid.HIPSPARSE_STATUS_ALLOC_FAILED – additional buffer for long rows could not be allocated.
HIPSPARSE_STATUS_NOT_SUPPORTED –
opA!= HIPSPARSE_OPERATION_NON_TRANSPOSE oropB!= HIPSPARSE_OPERATION_NON_TRANSPOSE.
hipsparseSpGEMMreuse_nnz()#
-
hipsparseStatus_t hipsparseSpGEMMreuse_nnz(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, hipsparseConstSpMatDescr_t matA, hipsparseConstSpMatDescr_t matB, hipsparseSpMatDescr_t matC, hipsparseSpGEMMAlg_t alg, hipsparseSpGEMMDescr_t spgemmDescr, size_t *bufferSize2, void *externalBuffer2, size_t *bufferSize3, void *externalBuffer3, size_t *bufferSize4, void *externalBuffer4)#
Nnz calculation step of the sparse matrix sparse matrix product:
\[ C' := \alpha \cdot op(A) \cdot op(B) + \beta \cdot C, \]where \(C'\), \(A\), \(B\), \(C\) are sparse matrices and \(C'\) and \(C\) have the same sparsity pattern.- Example (See full example below)
// Determine size of second, third, and fourth user allocated buffer hipsparseSpGEMMreuse_nnz(handle, opA, opB, matA, matB, matC, HIPSPARSE_SPGEMM_DEFAULT, spgemmDesc, &bufferSize2, NULL, &bufferSize3, NULL, &bufferSize4, NULL); hipMalloc((void**) &dBuffer2, bufferSize2); hipMalloc((void**) &dBuffer3, bufferSize3); hipMalloc((void**) &dBuffer4, bufferSize4); // COmpute sparsity pattern of C matrix and store in temporary buffers hipsparseSpGEMMreuse_nnz(handle, opA, opB, matA, matB, matC, HIPSPARSE_SPGEMM_DEFAULT, spgemmDesc, &bufferSize2, dBuffer2, &bufferSize3, dBuffer3, &bufferSize4, dBuffer4); // We can now free buffer 1 and 2 hipFree(dBuffer1); hipFree(dBuffer2);
- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] sparse matrix \(A\) operation type.
opB – [in] sparse matrix \(B\) operation type.
matA – [in] sparse matrix \(A\) descriptor.
matB – [in] sparse matrix \(B\) descriptor.
matC – [out] sparse matrix \(C\) descriptor.
alg – [in] SpGEMM algorithm for the SpGEMM computation.
spgemmDescr – [in] SpGEMM descriptor.
bufferSize2 – [out] number of bytes of the temporary storage
externalBuffer2.externalBuffer2 – [in] temporary storage buffer allocated by the user.
bufferSize3 – [out] number of bytes of the temporary storage
externalBuffer3.externalBuffer3 – [in] temporary storage buffer allocated by the user.
bufferSize4 – [out] number of bytes of the temporary storage
externalBuffer4.externalBuffer4 – [in] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,matA,matB,matC,bufferSize2,bufferSize3orbufferSize4pointer is invalid.HIPSPARSE_STATUS_ALLOC_FAILED – additional buffer for long rows could not be allocated.
HIPSPARSE_STATUS_NOT_SUPPORTED –
opA!= HIPSPARSE_OPERATION_NON_TRANSPOSE oropB!= HIPSPARSE_OPERATION_NON_TRANSPOSE.
hipsparseSpGEMMreuse_copy()#
-
hipsparseStatus_t hipsparseSpGEMMreuse_copy(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, hipsparseConstSpMatDescr_t matA, hipsparseConstSpMatDescr_t matB, hipsparseSpMatDescr_t matC, hipsparseSpGEMMAlg_t alg, hipsparseSpGEMMDescr_t spgemmDescr, size_t *bufferSize5, void *externalBuffer5)#
Copy step of the sparse matrix sparse matrix product:
\[ C' := \alpha \cdot op(A) \cdot op(B) + \beta \cdot C, \]where \(C'\), \(A\), \(B\), \(C\) are sparse matrices and \(C'\) and \(C\) have the same sparsity pattern.- Example (See full example below)
// Get matrix C non-zero entries nnzC int64_t rowsC, colsC, nnzC; hipsparseSpMatGetSize(matC, &rowsC, &colsC, &nnzC); // Allocate matrix C hipMalloc((void**) &dcsrColIndC, sizeof(int) * nnzC); hipMalloc((void**) &dcsrValC, sizeof(float) * nnzC); // Update matC with the new pointers. The C values array can be filled with data here // which is used if beta != 0. hipsparseCsrSetPointers(matC, dcsrRowPtrC, dcsrColIndC, dcsrValC); // Determine size of fifth user allocated buffer hipsparseSpGEMMreuse_copy(handle, opA, opB, matA, matB, matC, HIPSPARSE_SPGEMM_DEFAULT, spgemmDesc, &bufferSize5, NULL); hipMalloc((void**) &dBuffer5, bufferSize5); // Copy data from temporary buffers to the newly allocated C matrix hipsparseSpGEMMreuse_copy(handle, opA, opB, matA, matB, matC, HIPSPARSE_SPGEMM_DEFAULT, spgemmDesc, &bufferSize5, dBuffer5); // We can now free buffer 3 hipFree(dBuffer3);
- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] sparse matrix \(A\) operation type.
opB – [in] sparse matrix \(B\) operation type.
matA – [in] sparse matrix \(A\) descriptor.
matB – [in] sparse matrix \(B\) descriptor.
matC – [out] sparse matrix \(C\) descriptor.
alg – [in] SpGEMM algorithm for the SpGEMM computation.
spgemmDescr – [in] SpGEMM descriptor.
bufferSize5 – [out] number of bytes of the temporary storage
externalBuffer5.externalBuffer5 – [in] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,matA,matB,matC, orbufferSize5pointer is invalid.HIPSPARSE_STATUS_ALLOC_FAILED – additional buffer for long rows could not be allocated.
HIPSPARSE_STATUS_NOT_SUPPORTED –
opA!= HIPSPARSE_OPERATION_NON_TRANSPOSE oropB!= HIPSPARSE_OPERATION_NON_TRANSPOSE.
hipsparseSpGEMMreuse_compute()#
-
hipsparseStatus_t hipsparseSpGEMMreuse_compute(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstSpMatDescr_t matB, const void *beta, hipsparseSpMatDescr_t matC, hipDataType computeType, hipsparseSpGEMMAlg_t alg, hipsparseSpGEMMDescr_t spgemmDescr)#
Copy step of the sparse matrix sparse matrix product:
\[ C' := \alpha \cdot op(A) \cdot op(B) + \beta \cdot C, \]where \(C'\), \(A\), \(B\), \(C\) are sparse matrices and \(C'\) and \(C\) have the same sparsity pattern.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] sparse matrix \(A\) operation type.
opB – [in] sparse matrix \(B\) operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] sparse matrix \(A\) descriptor.
matB – [in] sparse matrix \(B\) descriptor.
beta – [in] scalar \(\beta\).
matC – [out] sparse matrix \(C\) descriptor.
computeType – [in] floating point precision for the SpGEMM computation.
alg – [in] SpGEMM algorithm for the SpGEMM computation.
spgemmDescr – [in] SpGEMM descriptor.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,beta,matA,matB, ormatCpointer is invalid.HIPSPARSE_STATUS_ALLOC_FAILED – additional buffer for long rows could not be allocated.
HIPSPARSE_STATUS_NOT_SUPPORTED –
opA!= HIPSPARSE_OPERATION_NON_TRANSPOSE oropB!= HIPSPARSE_OPERATION_NON_TRANSPOSE.
1int main(int argc, char* argv[])
2{
3 int m = 2;
4 int k = 2;
5 int n = 3;
6 int nnzA = 4;
7 int nnzB = 4;
8
9 float alpha{1.0f};
10 float beta{0.0f};
11
12 hipsparseOperation_t opA = HIPSPARSE_OPERATION_NON_TRANSPOSE;
13 hipsparseOperation_t opB = HIPSPARSE_OPERATION_NON_TRANSPOSE;
14 hipDataType computeType = HIP_R_32F;
15
16 // A, B, and C are m×k, k×n, and m×n
17
18 // A
19 std::vector<int> hcsrRowPtrA = {0, 2, 4};
20 std::vector<int> hcsrColIndA = {0, 1, 0, 1};
21 std::vector<float> hcsrValA = {1.0f, 2.0f, 3.0f, 4.0f};
22
23 // B
24 std::vector<int> hcsrRowPtrB = {0, 2, 4};
25 std::vector<int> hcsrColIndB = {1, 2, 0, 2};
26 std::vector<float> hcsrValB = {5.0f, 6.0f, 7.0f, 8.0f};
27
28 // Device memory management: Allocate and copy A, B
29 int* dcsrRowPtrA;
30 int* dcsrColIndA;
31 float* dcsrValA;
32 int* dcsrRowPtrB;
33 int* dcsrColIndB;
34 float* dcsrValB;
35 int* dcsrRowPtrC;
36 int* dcsrColIndC;
37 float* dcsrValC;
38 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrA, (m + 1) * sizeof(int)));
39 HIP_CHECK(hipMalloc((void**)&dcsrColIndA, nnzA * sizeof(int)));
40 HIP_CHECK(hipMalloc((void**)&dcsrValA, nnzA * sizeof(float)));
41 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrB, (k + 1) * sizeof(int)));
42 HIP_CHECK(hipMalloc((void**)&dcsrColIndB, nnzB * sizeof(int)));
43 HIP_CHECK(hipMalloc((void**)&dcsrValB, nnzB * sizeof(float)));
44 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrC, (m + 1) * sizeof(int)));
45
46 HIP_CHECK(
47 hipMemcpy(dcsrRowPtrA, hcsrRowPtrA.data(), (m + 1) * sizeof(int), hipMemcpyHostToDevice));
48 HIP_CHECK(
49 hipMemcpy(dcsrColIndA, hcsrColIndA.data(), nnzA * sizeof(int), hipMemcpyHostToDevice));
50 HIP_CHECK(hipMemcpy(dcsrValA, hcsrValA.data(), nnzA * sizeof(float), hipMemcpyHostToDevice));
51 HIP_CHECK(
52 hipMemcpy(dcsrRowPtrB, hcsrRowPtrB.data(), (k + 1) * sizeof(int), hipMemcpyHostToDevice));
53 HIP_CHECK(
54 hipMemcpy(dcsrColIndB, hcsrColIndB.data(), nnzB * sizeof(int), hipMemcpyHostToDevice));
55 HIP_CHECK(hipMemcpy(dcsrValB, hcsrValB.data(), nnzB * sizeof(float), hipMemcpyHostToDevice));
56
57 hipsparseHandle_t handle = NULL;
58 hipsparseSpMatDescr_t matA, matB, matC;
59 void* dBuffer1 = NULL;
60 void* dBuffer2 = NULL;
61 void* dBuffer3 = NULL;
62 void* dBuffer4 = NULL;
63 void* dBuffer5 = NULL;
64 size_t bufferSize1 = 0;
65 size_t bufferSize2 = 0;
66 size_t bufferSize3 = 0;
67 size_t bufferSize4 = 0;
68 size_t bufferSize5 = 0;
69
70 HIPSPARSE_CHECK(hipsparseCreate(&handle));
71
72 // Create sparse matrix A in CSR format
73 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
74 m,
75 k,
76 nnzA,
77 dcsrRowPtrA,
78 dcsrColIndA,
79 dcsrValA,
80 HIPSPARSE_INDEX_32I,
81 HIPSPARSE_INDEX_32I,
82 HIPSPARSE_INDEX_BASE_ZERO,
83 HIP_R_32F));
84 HIPSPARSE_CHECK(hipsparseCreateCsr(&matB,
85 k,
86 n,
87 nnzB,
88 dcsrRowPtrB,
89 dcsrColIndB,
90 dcsrValB,
91 HIPSPARSE_INDEX_32I,
92 HIPSPARSE_INDEX_32I,
93 HIPSPARSE_INDEX_BASE_ZERO,
94 HIP_R_32F));
95 HIPSPARSE_CHECK(hipsparseCreateCsr(&matC,
96 m,
97 n,
98 0,
99 dcsrRowPtrC,
100 NULL,
101 NULL,
102 HIPSPARSE_INDEX_32I,
103 HIPSPARSE_INDEX_32I,
104 HIPSPARSE_INDEX_BASE_ZERO,
105 HIP_R_32F));
106
107 hipsparseSpGEMMDescr_t spgemmDesc;
108 HIPSPARSE_CHECK(hipsparseSpGEMM_createDescr(&spgemmDesc));
109
110 // Determine size of first user allocated buffer
111 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_workEstimation(handle,
112 opA,
113 opB,
114 matA,
115 matB,
116 matC,
117 HIPSPARSE_SPGEMM_DEFAULT,
118 spgemmDesc,
119 &bufferSize1,
120 NULL));
121
122 HIP_CHECK(hipMalloc((void**)&dBuffer1, bufferSize1));
123
124 // Inspect the matrices A and B to determine the number of intermediate product in
125 // C = alpha * A * B
126 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_workEstimation(handle,
127 opA,
128 opB,
129 matA,
130 matB,
131 matC,
132 HIPSPARSE_SPGEMM_DEFAULT,
133 spgemmDesc,
134 &bufferSize1,
135 dBuffer1));
136
137 // Determine size of second, third, and fourth user allocated buffer
138 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_nnz(handle,
139 opA,
140 opB,
141 matA,
142 matB,
143 matC,
144 HIPSPARSE_SPGEMM_DEFAULT,
145 spgemmDesc,
146 &bufferSize2,
147 NULL,
148 &bufferSize3,
149 NULL,
150 &bufferSize4,
151 NULL));
152
153 HIP_CHECK(hipMalloc((void**)&dBuffer2, bufferSize2));
154 HIP_CHECK(hipMalloc((void**)&dBuffer3, bufferSize3));
155 HIP_CHECK(hipMalloc((void**)&dBuffer4, bufferSize4));
156
157 // Compute sparsity pattern of C matrix and store in temporary buffers
158 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_nnz(handle,
159 opA,
160 opB,
161 matA,
162 matB,
163 matC,
164 HIPSPARSE_SPGEMM_DEFAULT,
165 spgemmDesc,
166 &bufferSize2,
167 dBuffer2,
168 &bufferSize3,
169 dBuffer3,
170 &bufferSize4,
171 dBuffer4));
172
173 // We can now free buffer 1 and 2
174 HIP_CHECK(hipFree(dBuffer1));
175 HIP_CHECK(hipFree(dBuffer2));
176
177 // Get matrix C non-zero entries nnzC
178 int64_t rowsC, colsC, nnzC;
179 HIPSPARSE_CHECK(hipsparseSpMatGetSize(matC, &rowsC, &colsC, &nnzC));
180
181 // Allocate matrix C
182 HIP_CHECK(hipMalloc((void**)&dcsrColIndC, sizeof(int) * nnzC));
183 HIP_CHECK(hipMalloc((void**)&dcsrValC, sizeof(float) * nnzC));
184
185 // Update matC with the new pointers. The C values array can be filled with data here
186 // which is used if beta != 0.
187 HIPSPARSE_CHECK(hipsparseCsrSetPointers(matC, dcsrRowPtrC, dcsrColIndC, dcsrValC));
188
189 // Determine size of fifth user allocated buffer
190 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_copy(handle,
191 opA,
192 opB,
193 matA,
194 matB,
195 matC,
196 HIPSPARSE_SPGEMM_DEFAULT,
197 spgemmDesc,
198 &bufferSize5,
199 NULL));
200
201 HIP_CHECK(hipMalloc((void**)&dBuffer5, bufferSize5));
202
203 // Copy data from temporary buffers to the newly allocated C matrix
204 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_copy(handle,
205 opA,
206 opB,
207 matA,
208 matB,
209 matC,
210 HIPSPARSE_SPGEMM_DEFAULT,
211 spgemmDesc,
212 &bufferSize5,
213 dBuffer5));
214
215 // We can now free buffer 3
216 HIP_CHECK(hipFree(dBuffer3));
217
218 // Compute C' = alpha * A * B + beta * C
219 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_compute(handle,
220 opA,
221 opB,
222 &alpha,
223 matA,
224 matB,
225 &beta,
226 matC,
227 computeType,
228 HIPSPARSE_SPGEMM_DEFAULT,
229 spgemmDesc));
230
231 // Copy results back to host if required
232 std::vector<int> hcsrRowPtrC(m + 1);
233 std::vector<int> hcsrColIndC(nnzC);
234 std::vector<float> hcsrValC(nnzC);
235 HIP_CHECK(
236 hipMemcpy(hcsrRowPtrC.data(), dcsrRowPtrC, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
237 HIP_CHECK(
238 hipMemcpy(hcsrColIndC.data(), dcsrColIndC, sizeof(int) * nnzC, hipMemcpyDeviceToHost));
239 HIP_CHECK(hipMemcpy(hcsrValC.data(), dcsrValC, sizeof(float) * nnzC, hipMemcpyDeviceToHost));
240
241 // Update dcsrValA, dcsrValB with new values
242 for(size_t i = 0; i < hcsrValA.size(); i++)
243 {
244 hcsrValA[i] = 1.0f;
245 }
246 for(size_t i = 0; i < hcsrValB.size(); i++)
247 {
248 hcsrValB[i] = 2.0f;
249 }
250
251 HIP_CHECK(hipMemcpy(dcsrValA, hcsrValA.data(), sizeof(float) * nnzA, hipMemcpyHostToDevice));
252 HIP_CHECK(hipMemcpy(dcsrValB, hcsrValB.data(), sizeof(float) * nnzB, hipMemcpyHostToDevice));
253
254 // Compute C' = alpha * A * B + beta * C again with the new A and B values
255 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_compute(handle,
256 opA,
257 opB,
258 &alpha,
259 matA,
260 matB,
261 &beta,
262 matC,
263 computeType,
264 HIPSPARSE_SPGEMM_DEFAULT,
265 spgemmDesc));
266
267 // Copy results back to host if required
268 HIP_CHECK(
269 hipMemcpy(hcsrRowPtrC.data(), dcsrRowPtrC, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
270 HIP_CHECK(
271 hipMemcpy(hcsrColIndC.data(), dcsrColIndC, sizeof(int) * nnzC, hipMemcpyDeviceToHost));
272 HIP_CHECK(hipMemcpy(hcsrValC.data(), dcsrValC, sizeof(float) * nnzC, hipMemcpyDeviceToHost));
273
274 // Destroy matrix descriptors and handles
275 HIPSPARSE_CHECK(hipsparseSpGEMM_destroyDescr(spgemmDesc));
276 HIPSPARSE_CHECK(hipsparseDestroySpMat(matA));
277 HIPSPARSE_CHECK(hipsparseDestroySpMat(matB));
278 HIPSPARSE_CHECK(hipsparseDestroySpMat(matC));
279 HIPSPARSE_CHECK(hipsparseDestroy(handle));
280
281 // Free device memory
282 HIP_CHECK(hipFree(dBuffer4));
283 HIP_CHECK(hipFree(dBuffer5));
284 HIP_CHECK(hipFree(dcsrRowPtrA));
285 HIP_CHECK(hipFree(dcsrColIndA));
286 HIP_CHECK(hipFree(dcsrValA));
287 HIP_CHECK(hipFree(dcsrRowPtrB));
288 HIP_CHECK(hipFree(dcsrColIndB));
289 HIP_CHECK(hipFree(dcsrValB));
290 HIP_CHECK(hipFree(dcsrRowPtrC));
291 HIP_CHECK(hipFree(dcsrColIndC));
292 HIP_CHECK(hipFree(dcsrValC));
293
294 return 0;
295}
1int main(int argc, char* argv[])
2{
3 int m = 2;
4 int k = 2;
5 int n = 3;
6 int nnzA = 4;
7 int nnzB = 4;
8
9 float alpha = 1.0;
10 float beta = 0.0;
11
12 hipsparseOperation_t opA = HIPSPARSE_OPERATION_NON_TRANSPOSE;
13 hipsparseOperation_t opB = HIPSPARSE_OPERATION_NON_TRANSPOSE;
14 hipDataType computeType = HIP_R_32F;
15
16 // A, B, and C are m×k, k×n, and m×n
17
18 // A
19 int hcsrRowPtrA[] = {0, 2, 4};
20 int hcsrColIndA[] = {0, 1, 0, 1};
21 float hcsrValA[] = {1.0, 2.0, 3.0, 4.0};
22
23 // B
24 int hcsrRowPtrB[] = {0, 2, 4};
25 int hcsrColIndB[] = {1, 2, 0, 2};
26 float hcsrValB[] = {5.0, 6.0, 7.0, 8.0};
27
28 // Device memory management: Allocate and copy A, B
29 int* dcsrRowPtrA;
30 int* dcsrColIndA;
31 float* dcsrValA;
32 int* dcsrRowPtrB;
33 int* dcsrColIndB;
34 float* dcsrValB;
35 int* dcsrRowPtrC;
36 int* dcsrColIndC;
37 float* dcsrValC;
38 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrA, (m + 1) * sizeof(int)));
39 HIP_CHECK(hipMalloc((void**)&dcsrColIndA, nnzA * sizeof(int)));
40 HIP_CHECK(hipMalloc((void**)&dcsrValA, nnzA * sizeof(float)));
41 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrB, (k + 1) * sizeof(int)));
42 HIP_CHECK(hipMalloc((void**)&dcsrColIndB, nnzB * sizeof(int)));
43 HIP_CHECK(hipMalloc((void**)&dcsrValB, nnzB * sizeof(float)));
44 HIP_CHECK(hipMalloc((void**)&dcsrRowPtrC, (m + 1) * sizeof(int)));
45
46 HIP_CHECK(hipMemcpy(dcsrRowPtrA, hcsrRowPtrA, (m + 1) * sizeof(int), hipMemcpyHostToDevice));
47 HIP_CHECK(hipMemcpy(dcsrColIndA, hcsrColIndA, nnzA * sizeof(int), hipMemcpyHostToDevice));
48 HIP_CHECK(hipMemcpy(dcsrValA, hcsrValA, nnzA * sizeof(float), hipMemcpyHostToDevice));
49 HIP_CHECK(hipMemcpy(dcsrRowPtrB, hcsrRowPtrB, (k + 1) * sizeof(int), hipMemcpyHostToDevice));
50 HIP_CHECK(hipMemcpy(dcsrColIndB, hcsrColIndB, nnzB * sizeof(int), hipMemcpyHostToDevice));
51 HIP_CHECK(hipMemcpy(dcsrValB, hcsrValB, nnzB * sizeof(float), hipMemcpyHostToDevice));
52
53 hipsparseHandle_t handle = NULL;
54 hipsparseSpMatDescr_t matA, matB, matC;
55 void* dBuffer1 = NULL;
56 void* dBuffer2 = NULL;
57 void* dBuffer3 = NULL;
58 void* dBuffer4 = NULL;
59 void* dBuffer5 = NULL;
60 size_t bufferSize1 = 0;
61 size_t bufferSize2 = 0;
62 size_t bufferSize3 = 0;
63 size_t bufferSize4 = 0;
64 size_t bufferSize5 = 0;
65
66 HIPSPARSE_CHECK(hipsparseCreate(&handle));
67
68 // Create sparse matrix A in CSR format
69 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
70 m,
71 k,
72 nnzA,
73 dcsrRowPtrA,
74 dcsrColIndA,
75 dcsrValA,
76 HIPSPARSE_INDEX_32I,
77 HIPSPARSE_INDEX_32I,
78 HIPSPARSE_INDEX_BASE_ZERO,
79 HIP_R_32F));
80 HIPSPARSE_CHECK(hipsparseCreateCsr(&matB,
81 k,
82 n,
83 nnzB,
84 dcsrRowPtrB,
85 dcsrColIndB,
86 dcsrValB,
87 HIPSPARSE_INDEX_32I,
88 HIPSPARSE_INDEX_32I,
89 HIPSPARSE_INDEX_BASE_ZERO,
90 HIP_R_32F));
91 HIPSPARSE_CHECK(hipsparseCreateCsr(&matC,
92 m,
93 n,
94 0,
95 dcsrRowPtrC,
96 NULL,
97 NULL,
98 HIPSPARSE_INDEX_32I,
99 HIPSPARSE_INDEX_32I,
100 HIPSPARSE_INDEX_BASE_ZERO,
101 HIP_R_32F));
102
103 hipsparseSpGEMMDescr_t spgemmDesc;
104 HIPSPARSE_CHECK(hipsparseSpGEMM_createDescr(&spgemmDesc));
105
106 // Determine size of first user allocated buffer
107 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_workEstimation(handle,
108 opA,
109 opB,
110 matA,
111 matB,
112 matC,
113 HIPSPARSE_SPGEMM_DEFAULT,
114 spgemmDesc,
115 &bufferSize1,
116 NULL));
117
118 HIP_CHECK(hipMalloc((void**)&dBuffer1, bufferSize1));
119
120 // Inspect the matrices A and B to determine the number of intermediate product in
121 // C = alpha * A * B
122 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_workEstimation(handle,
123 opA,
124 opB,
125 matA,
126 matB,
127 matC,
128 HIPSPARSE_SPGEMM_DEFAULT,
129 spgemmDesc,
130 &bufferSize1,
131 dBuffer1));
132
133 // Determine size of second, third, and fourth user allocated buffer
134 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_nnz(handle,
135 opA,
136 opB,
137 matA,
138 matB,
139 matC,
140 HIPSPARSE_SPGEMM_DEFAULT,
141 spgemmDesc,
142 &bufferSize2,
143 NULL,
144 &bufferSize3,
145 NULL,
146 &bufferSize4,
147 NULL));
148
149 HIP_CHECK(hipMalloc((void**)&dBuffer2, bufferSize2));
150 HIP_CHECK(hipMalloc((void**)&dBuffer3, bufferSize3));
151 HIP_CHECK(hipMalloc((void**)&dBuffer4, bufferSize4));
152
153 // Compute sparsity pattern of C matrix and store in temporary buffers
154 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_nnz(handle,
155 opA,
156 opB,
157 matA,
158 matB,
159 matC,
160 HIPSPARSE_SPGEMM_DEFAULT,
161 spgemmDesc,
162 &bufferSize2,
163 dBuffer2,
164 &bufferSize3,
165 dBuffer3,
166 &bufferSize4,
167 dBuffer4));
168
169 // We can now free buffer 1 and 2
170 HIP_CHECK(hipFree(dBuffer1));
171 HIP_CHECK(hipFree(dBuffer2));
172
173 // Get matrix C non-zero entries nnzC
174 int64_t rowsC, colsC, nnzC;
175 HIPSPARSE_CHECK(hipsparseSpMatGetSize(matC, &rowsC, &colsC, &nnzC));
176
177 // Allocate matrix C
178 HIP_CHECK(hipMalloc((void**)&dcsrColIndC, sizeof(int) * nnzC));
179 HIP_CHECK(hipMalloc((void**)&dcsrValC, sizeof(float) * nnzC));
180
181 // Update matC with the new pointers. The C values array can be filled with data here
182 // which is used if beta != 0.
183 HIPSPARSE_CHECK(hipsparseCsrSetPointers(matC, dcsrRowPtrC, dcsrColIndC, dcsrValC));
184
185 // Determine size of fifth user allocated buffer
186 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_copy(handle,
187 opA,
188 opB,
189 matA,
190 matB,
191 matC,
192 HIPSPARSE_SPGEMM_DEFAULT,
193 spgemmDesc,
194 &bufferSize5,
195 NULL));
196
197 HIP_CHECK(hipMalloc((void**)&dBuffer5, bufferSize5));
198
199 // Copy data from temporary buffers to the newly allocated C matrix
200 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_copy(handle,
201 opA,
202 opB,
203 matA,
204 matB,
205 matC,
206 HIPSPARSE_SPGEMM_DEFAULT,
207 spgemmDesc,
208 &bufferSize5,
209 dBuffer5));
210
211 // We can now free buffer 3
212 HIP_CHECK(hipFree(dBuffer3));
213
214 // Compute C' = alpha * A * B + beta * C
215 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_compute(handle,
216 opA,
217 opB,
218 &alpha,
219 matA,
220 matB,
221 &beta,
222 matC,
223 computeType,
224 HIPSPARSE_SPGEMM_DEFAULT,
225 spgemmDesc));
226
227 // Copy results back to host if required
228 int* hcsrRowPtrC = (int*)malloc((m + 1) * sizeof(int));
229 int* hcsrColIndC = (int*)malloc((nnzC) * sizeof(int));
230 float hcsrValC[nnzC];
231 HIP_CHECK(hipMemcpy(hcsrRowPtrC, dcsrRowPtrC, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
232 HIP_CHECK(hipMemcpy(hcsrColIndC, dcsrColIndC, sizeof(int) * nnzC, hipMemcpyDeviceToHost));
233 HIP_CHECK(hipMemcpy(hcsrValC, dcsrValC, sizeof(float) * nnzC, hipMemcpyDeviceToHost));
234
235 // Update dcsrValA, dcsrValB with new values
236 for(size_t i = 0; i < nnzA; i++)
237 {
238 hcsrValA[i] = 1.0;
239 }
240 for(size_t i = 0; i < nnzB; i++)
241 {
242 hcsrValB[i] = 2.0;
243 }
244
245 HIP_CHECK(hipMemcpy(dcsrValA, hcsrValA, sizeof(float) * nnzA, hipMemcpyHostToDevice));
246 HIP_CHECK(hipMemcpy(dcsrValB, hcsrValB, sizeof(float) * nnzB, hipMemcpyHostToDevice));
247
248 // Compute C' = alpha * A * B + beta * C again with the new A and B values
249 HIPSPARSE_CHECK(hipsparseSpGEMMreuse_compute(handle,
250 opA,
251 opB,
252 &alpha,
253 matA,
254 matB,
255 &beta,
256 matC,
257 computeType,
258 HIPSPARSE_SPGEMM_DEFAULT,
259 spgemmDesc));
260
261 // Copy results back to host if required
262 HIP_CHECK(hipMemcpy(hcsrRowPtrC, dcsrRowPtrC, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
263 HIP_CHECK(hipMemcpy(hcsrColIndC, dcsrColIndC, sizeof(int) * nnzC, hipMemcpyDeviceToHost));
264 HIP_CHECK(hipMemcpy(hcsrValC, dcsrValC, sizeof(float) * nnzC, hipMemcpyDeviceToHost));
265
266 // Destroy matrix descriptors and handles
267 HIPSPARSE_CHECK(hipsparseSpGEMM_destroyDescr(spgemmDesc));
268 HIPSPARSE_CHECK(hipsparseDestroySpMat(matA));
269 HIPSPARSE_CHECK(hipsparseDestroySpMat(matB));
270 HIPSPARSE_CHECK(hipsparseDestroySpMat(matC));
271 HIPSPARSE_CHECK(hipsparseDestroy(handle));
272
273 // Free device memory
274 HIP_CHECK(hipFree(dBuffer4));
275 HIP_CHECK(hipFree(dBuffer5));
276 HIP_CHECK(hipFree(dcsrRowPtrA));
277 HIP_CHECK(hipFree(dcsrColIndA));
278 HIP_CHECK(hipFree(dcsrValA));
279 HIP_CHECK(hipFree(dcsrRowPtrB));
280 HIP_CHECK(hipFree(dcsrColIndB));
281 HIP_CHECK(hipFree(dcsrValB));
282 HIP_CHECK(hipFree(dcsrRowPtrC));
283 HIP_CHECK(hipFree(dcsrColIndC));
284 HIP_CHECK(hipFree(dcsrValC));
285
286 return 0;
287}
hipsparseSDDMM_bufferSize()#
-
hipsparseStatus_t hipsparseSDDMM_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstDnMatDescr_t A, hipsparseConstDnMatDescr_t B, const void *beta, hipsparseSpMatDescr_t C, hipDataType computeType, hipsparseSDDMMAlg_t alg, size_t *pBufferSizeInBytes)#
hipsparseSDDMM_bufferSizereturns the size of the required buffer needed when computing the sampled dense-dense matrix multiplication:\[ C := \alpha (op(A) \cdot op(B)) \circ spy(C) + \beta \cdot C, \]where \(C\) is a sparse matrix and \(A\) and \(B\) are dense matrices. This routine is used in conjunction with hipsparseSDDMM_preprocess() and hipsparseSDDMM().hipsparseSDDMM_bufferSizesupports multiple combinations of data types and compute types. See hipsparseSDDMM for a complete listing of all the data type and compute type combinations available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] dense matrix \(A\) operation type.
opB – [in] dense matrix \(B\) operation type.
alpha – [in] scalar \(\alpha\).
A – [in] dense matrix \(A\) descriptor.
B – [in] dense matrix \(B\) descriptor.
beta – [in] scalar \(\beta\).
C – [inout] sparse matrix \(C\) descriptor.
computeType – [in] floating point precision for the SDDMM computation.
alg – [in] specification of the algorithm to use.
pBufferSizeInBytes – [out] number of bytes of the temporary storage buffer.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_NOT_INITIALIZED –
handleis not initialized.HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,beta,A,B,CorpBufferSizeInBytesis nullptr, oropAoropBis invalid.HIPSPARSE_STATUS_NOT_SUPPORTED –
opA== HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE oropB== HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE.
hipsparseSDDMM_preprocess()#
-
hipsparseStatus_t hipsparseSDDMM_preprocess(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstDnMatDescr_t A, hipsparseConstDnMatDescr_t B, const void *beta, hipsparseSpMatDescr_t C, hipDataType computeType, hipsparseSDDMMAlg_t alg, void *tempBuffer)#
hipsparseSDDMM_preprocessperforms the required preprocessing used when computing the sampled dense dense matrix multiplication:\[ C := \alpha (op(A) \cdot op(B)) \circ spy(C) + \beta \cdot C, \]where \(C\) is a sparse matrix and \(A\) and \(B\) are dense matrices. This routine is used in conjunction with hipsparseSDDMM().hipsparseSDDMM_preprocesssupports multiple combinations of data types and compute types. See hipsparseSDDMM for a complete listing of all the data type and compute type combinations available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] dense matrix \(A\) operation type.
opB – [in] dense matrix \(B\) operation type.
alpha – [in] scalar \(\alpha\).
A – [in] dense matrix \(A\) descriptor.
B – [in] dense matrix \(B\) descriptor.
beta – [in] scalar \(\beta\).
C – [inout] sparse matrix \(C\) descriptor.
computeType – [in] floating point precision for the SDDMM computation.
alg – [in] specification of the algorithm to use.
tempBuffer – [in] temporary storage buffer allocated by the user. The size must be greater or equal to the size obtained with hipsparseSDDMM_bufferSize.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,beta,A,B,CortempBufferpointer is invalid or the value ofopAoropBis incorrect.HIPSPARSE_STATUS_NOT_SUPPORTED –
opA== HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE oropB== HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE.
hipsparseSDDMM()#
-
hipsparseStatus_t hipsparseSDDMM(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstDnMatDescr_t A, hipsparseConstDnMatDescr_t B, const void *beta, hipsparseSpMatDescr_t C, hipDataType computeType, hipsparseSDDMMAlg_t alg, void *tempBuffer)#
Sampled Dense-Dense Matrix Multiplication.
hipsparseSDDMMmultiplies the scalar \(\alpha\) with the dense \(m \times k\) matrix \(op(A)\), the dense \(k \times n\) matrix \(op(B)\), filtered by the sparsity pattern of the \(m \times n\) sparse matrix \(C\) and adds the result to \(C\) scaled by \(\beta\). The final result is stored in the sparse \(m \times n\) matrix \(C\), such that\[ C := \alpha ( op(A) \cdot op(B) ) \circ spy(C) + \beta C, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if op(A) == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ A^T, & \text{if op(A) == HIPSPARSE_OPERATION_TRANSPOSE} \\ \end{array} \right. \end{split}\],\[\begin{split} op(B) = \left\{ \begin{array}{ll} B, & \text{if op(B) == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ B^T, & \text{if op(B) == HIPSPARSE_OPERATION_TRANSPOSE} \\ \end{array} \right. \end{split}\]and\[\begin{split} spy(C)_{ij} = \left\{ \begin{array}{ll} 1, & \text{ if C_{ij} != 0} \\ 0, & \text{ otherwise} \\ \end{array} \right. \end{split}\]Computing the above sampled dense-dense multiplication requires three steps to complete. First, the user calls hipsparseSDDMM_bufferSize to determine the size of the required temporary storage buffer. Next, the user allocates this buffer and calls hipsparseSDDMM_preprocess which performs any analysis of the input matrices that may be required. Finally, the user calls
hipsparseSDDMMto complete the computation. Once all calls tohipsparseSDDMMare complete, the temporary buffer can be deallocated.hipsparseSDDMMsupports different algorithms which can provide better performance for different matrices.CSR/CSC Algorithms
HIPSPARSE_SDDMM_ALG_DEFAULT
Currently,
hipsparseSDDMMonly supports the uniform precisions indicated in the table below. For the sparse matrix \(C\),hipsparseSDDMMsupports the index types HIPSPARSE_INDEX_32I and HIPSPARSE_INDEX_64I.- Uniform Precisions:
A / B / C / compute_type
HIP_R_16F
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
- Mixed precisions:
A / B
C
compute_type
HIP_R_16F
HIP_R_32F
HIP_R_32F
HIP_R_16F
HIP_R_16F
HIP_R_32F
- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] dense matrix \(A\) operation type.
opB – [in] dense matrix \(B\) operation type.
alpha – [in] scalar \(\alpha\).
A – [in] dense matrix \(A\) descriptor.
B – [in] dense matrix \(B\) descriptor.
beta – [in] scalar \(\beta\).
C – [inout] sparse matrix \(C\) descriptor.
computeType – [in] floating point precision for the SDDMM computation.
alg – [in] specification of the algorithm to use.
tempBuffer – [in] temporary storage buffer allocated by the user. The size must be greater or equal to the size obtained with hipsparseSDDMM_bufferSize.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,beta,A,B,CortempBufferpointer is invalid or the value ofopAoropBis incorrect.HIPSPARSE_STATUS_NOT_SUPPORTED –
opA== HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE oropB== HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE.
1int main(int argc, char* argv[])
2{
3 // hipSPARSE handle
4 hipsparseHandle_t handle;
5 HIPSPARSE_CHECK(hipsparseCreate(&handle));
6
7 __half halpha = 1.0;
8 __half hbeta = 0.0;
9
10 // A, B, and C are mxk, kxn, and mxn
11 int m = 4;
12 int k = 3;
13 int n = 2;
14 int nnzC = 5;
15
16 // 2 3 -1
17 // A = 0 2 1
18 // 0 0 5
19 // 0 -2 0.5
20
21 // 0 4
22 // B = 1 0
23 // -2 0.5
24
25 // 1 0 1 0
26 // C = 2 3 spy(C) = 1 1
27 // 0 0 0 0
28 // 4 5 1 1
29
30 std::vector<__half> hA = {2.0, 3.0, -1.0, 0.0, 2.0, 1.0, 0.0, 0.0, 5.0, 0.0, -2.0, 0.5};
31 std::vector<__half> hB = {0.0, 4.0, 1.0, 0.0, -2.0, 0.5};
32
33 std::vector<int> hcsr_row_ptrC = {0, 1, 3, 3, 5};
34 std::vector<int> hcsr_col_indC = {0, 0, 1, 0, 1};
35 std::vector<__half> hcsr_valC = {1.0, 2.0, 3.0, 4.0, 5.0};
36
37 __half* dA = nullptr;
38 __half* dB = nullptr;
39 HIP_CHECK(hipMalloc((void**)&dA, sizeof(__half) * m * k));
40 HIP_CHECK(hipMalloc((void**)&dB, sizeof(__half) * k * n));
41
42 int* dcsr_row_ptrC = nullptr;
43 int* dcsr_col_indC = nullptr;
44 __half* dcsr_valC = nullptr;
45 HIP_CHECK(hipMalloc((void**)&dcsr_row_ptrC, sizeof(int) * (m + 1)));
46 HIP_CHECK(hipMalloc((void**)&dcsr_col_indC, sizeof(int) * nnzC));
47 HIP_CHECK(hipMalloc((void**)&dcsr_valC, sizeof(__half) * nnzC));
48
49 HIP_CHECK(hipMemcpy(dA, hA.data(), sizeof(__half) * m * k, hipMemcpyHostToDevice));
50 HIP_CHECK(hipMemcpy(dB, hB.data(), sizeof(__half) * k * n, hipMemcpyHostToDevice));
51
52 HIP_CHECK(hipMemcpy(
53 dcsr_row_ptrC, hcsr_row_ptrC.data(), sizeof(int) * (m + 1), hipMemcpyHostToDevice));
54 HIP_CHECK(
55 hipMemcpy(dcsr_col_indC, hcsr_col_indC.data(), sizeof(int) * nnzC, hipMemcpyHostToDevice));
56 HIP_CHECK(hipMemcpy(dcsr_valC, hcsr_valC.data(), sizeof(__half) * nnzC, hipMemcpyHostToDevice));
57
58 hipsparseDnMatDescr_t matA;
59 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matA, m, k, k, dA, HIP_R_16F, HIPSPARSE_ORDER_ROW));
60
61 hipsparseDnMatDescr_t matB;
62 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matB, k, n, n, dB, HIP_R_16F, HIPSPARSE_ORDER_ROW));
63
64 hipsparseSpMatDescr_t matC;
65 HIPSPARSE_CHECK(hipsparseCreateCsr(&matC,
66 m,
67 n,
68 nnzC,
69 dcsr_row_ptrC,
70 dcsr_col_indC,
71 dcsr_valC,
72 HIPSPARSE_INDEX_32I,
73 HIPSPARSE_INDEX_32I,
74 HIPSPARSE_INDEX_BASE_ZERO,
75 HIP_R_16F));
76
77 size_t buffer_size = 0;
78 HIPSPARSE_CHECK(hipsparseSDDMM_bufferSize(handle,
79 HIPSPARSE_OPERATION_NON_TRANSPOSE,
80 HIPSPARSE_OPERATION_NON_TRANSPOSE,
81 &halpha,
82 matA,
83 matB,
84 &hbeta,
85 matC,
86 HIP_R_16F,
87 HIPSPARSE_SDDMM_ALG_DEFAULT,
88 &buffer_size));
89
90 void* dbuffer = nullptr;
91 HIP_CHECK(hipMalloc((void**)&dbuffer, buffer_size));
92
93 HIPSPARSE_CHECK(hipsparseSDDMM_preprocess(handle,
94 HIPSPARSE_OPERATION_NON_TRANSPOSE,
95 HIPSPARSE_OPERATION_NON_TRANSPOSE,
96 &halpha,
97 matA,
98 matB,
99 &hbeta,
100 matC,
101 HIP_R_16F,
102 HIPSPARSE_SDDMM_ALG_DEFAULT,
103 dbuffer));
104
105 HIPSPARSE_CHECK(hipsparseSDDMM(handle,
106 HIPSPARSE_OPERATION_NON_TRANSPOSE,
107 HIPSPARSE_OPERATION_NON_TRANSPOSE,
108 &halpha,
109 matA,
110 matB,
111 &hbeta,
112 matC,
113 HIP_R_16F,
114 HIPSPARSE_SDDMM_ALG_DEFAULT,
115 dbuffer));
116
117 HIP_CHECK(hipMemcpy(
118 hcsr_row_ptrC.data(), dcsr_row_ptrC, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
119 HIP_CHECK(
120 hipMemcpy(hcsr_col_indC.data(), dcsr_col_indC, sizeof(int) * nnzC, hipMemcpyDeviceToHost));
121 HIP_CHECK(hipMemcpy(hcsr_valC.data(), dcsr_valC, sizeof(__half) * nnzC, hipMemcpyDeviceToHost));
122
123 std::cout << "C" << std::endl;
124 for(int i = 0; i < m; i++)
125 {
126 int start = hcsr_row_ptrC[i];
127 int end = hcsr_row_ptrC[i + 1];
128
129 std::vector<__half> temp(n, 0.0);
130 for(int j = start; j < end; j++)
131 {
132 temp[hcsr_col_indC[j]] = hcsr_valC[j];
133 }
134
135 for(int j = 0; j < n; j++)
136 {
137 std::cout << static_cast<float>(temp[j]) << " ";
138 }
139 std::cout << std::endl;
140 }
141 std::cout << std::endl;
142
143 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matA));
144 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matB));
145 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matC));
146 HIPSPARSE_CHECK(hipsparseDestroy(handle));
147
148 HIP_CHECK(hipFree(dA));
149 HIP_CHECK(hipFree(dB));
150 HIP_CHECK(hipFree(dcsr_row_ptrC));
151 HIP_CHECK(hipFree(dcsr_col_indC));
152 HIP_CHECK(hipFree(dcsr_valC));
153 HIP_CHECK(hipFree(dbuffer));
154
155 return 0;
156}
1int main(int argc, char* argv[])
2{
3 // hipSPARSE handle
4 hipsparseHandle_t handle;
5 HIPSPARSE_CHECK(hipsparseCreate(&handle));
6
7 float halpha = 1.0;
8 float hbeta = 0.0;
9
10 // A, B, and C are mxk, kxn, and mxn
11 int m = 4;
12 int k = 3;
13 int n = 2;
14 int nnzC = 5;
15
16 // 2 3 -1
17 // A = 0 2 1
18 // 0 0 5
19 // 0 -2 0.5
20
21 // 0 4
22 // B = 1 0
23 // -2 0.5
24
25 // 1 0 1 0
26 // C = 2 3 spy(C) = 1 1
27 // 0 0 0 0
28 // 4 5 1 1
29
30 float hA[] = {2.0, 3.0, -1.0, 0.0, 2.0, 1.0, 0.0, 0.0, 5.0, 0.0, -2.0, 0.5};
31 float hB[] = {0.0, 4.0, 1.0, 0.0, -2.0, 0.5};
32
33 int hcsr_row_ptrC[] = {0, 1, 3, 3, 5};
34 int hcsr_col_indC[] = {0, 0, 1, 0, 1};
35 float hcsr_valC[] = {1.0, 2.0, 3.0, 4.0, 5.0};
36
37 float* dA = NULL;
38 float* dB = NULL;
39 HIP_CHECK(hipMalloc((void**)&dA, sizeof(float) * m * k));
40 HIP_CHECK(hipMalloc((void**)&dB, sizeof(float) * k * n));
41
42 int* dcsr_row_ptrC = NULL;
43 int* dcsr_col_indC = NULL;
44 float* dcsr_valC = NULL;
45 HIP_CHECK(hipMalloc((void**)&dcsr_row_ptrC, sizeof(int) * (m + 1)));
46 HIP_CHECK(hipMalloc((void**)&dcsr_col_indC, sizeof(int) * nnzC));
47 HIP_CHECK(hipMalloc((void**)&dcsr_valC, sizeof(float) * nnzC));
48
49 HIP_CHECK(hipMemcpy(dA, hA, sizeof(float) * m * k, hipMemcpyHostToDevice));
50 HIP_CHECK(hipMemcpy(dB, hB, sizeof(float) * k * n, hipMemcpyHostToDevice));
51
52 HIP_CHECK(
53 hipMemcpy(dcsr_row_ptrC, hcsr_row_ptrC, sizeof(int) * (m + 1), hipMemcpyHostToDevice));
54 HIP_CHECK(hipMemcpy(dcsr_col_indC, hcsr_col_indC, sizeof(int) * nnzC, hipMemcpyHostToDevice));
55 HIP_CHECK(hipMemcpy(dcsr_valC, hcsr_valC, sizeof(float) * nnzC, hipMemcpyHostToDevice));
56
57 hipsparseDnMatDescr_t matA;
58 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matA, m, k, k, dA, HIP_R_32F, HIPSPARSE_ORDER_ROW));
59
60 hipsparseDnMatDescr_t matB;
61 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matB, k, n, n, dB, HIP_R_32F, HIPSPARSE_ORDER_ROW));
62
63 hipsparseSpMatDescr_t matC;
64 HIPSPARSE_CHECK(hipsparseCreateCsr(&matC,
65 m,
66 n,
67 nnzC,
68 dcsr_row_ptrC,
69 dcsr_col_indC,
70 dcsr_valC,
71 HIPSPARSE_INDEX_32I,
72 HIPSPARSE_INDEX_32I,
73 HIPSPARSE_INDEX_BASE_ZERO,
74 HIP_R_32F));
75
76 size_t buffer_size = 0;
77 HIPSPARSE_CHECK(hipsparseSDDMM_bufferSize(handle,
78 HIPSPARSE_OPERATION_NON_TRANSPOSE,
79 HIPSPARSE_OPERATION_NON_TRANSPOSE,
80 &halpha,
81 matA,
82 matB,
83 &hbeta,
84 matC,
85 HIP_R_32F,
86 HIPSPARSE_SDDMM_ALG_DEFAULT,
87 &buffer_size));
88
89 void* dbuffer = NULL;
90 HIP_CHECK(hipMalloc((void**)&dbuffer, buffer_size));
91
92 HIPSPARSE_CHECK(hipsparseSDDMM_preprocess(handle,
93 HIPSPARSE_OPERATION_NON_TRANSPOSE,
94 HIPSPARSE_OPERATION_NON_TRANSPOSE,
95 &halpha,
96 matA,
97 matB,
98 &hbeta,
99 matC,
100 HIP_R_32F,
101 HIPSPARSE_SDDMM_ALG_DEFAULT,
102 dbuffer));
103
104 HIPSPARSE_CHECK(hipsparseSDDMM(handle,
105 HIPSPARSE_OPERATION_NON_TRANSPOSE,
106 HIPSPARSE_OPERATION_NON_TRANSPOSE,
107 &halpha,
108 matA,
109 matB,
110 &hbeta,
111 matC,
112 HIP_R_32F,
113 HIPSPARSE_SDDMM_ALG_DEFAULT,
114 dbuffer));
115
116 HIP_CHECK(
117 hipMemcpy(hcsr_row_ptrC, dcsr_row_ptrC, sizeof(int) * (m + 1), hipMemcpyDeviceToHost));
118 HIP_CHECK(hipMemcpy(hcsr_col_indC, dcsr_col_indC, sizeof(int) * nnzC, hipMemcpyDeviceToHost));
119 HIP_CHECK(hipMemcpy(hcsr_valC, dcsr_valC, sizeof(float) * nnzC, hipMemcpyDeviceToHost));
120
121 printf("C\n");
122 for(int i = 0; i < m; i++)
123 {
124 int start = hcsr_row_ptrC[i];
125 int end = hcsr_row_ptrC[i + 1];
126
127 float* temp = (float*)malloc(n * sizeof(float));
128 for(int j = start; j < end; j++)
129 {
130 temp[hcsr_col_indC[j]] = hcsr_valC[j];
131 }
132
133 for(int j = 0; j < n; j++)
134 {
135 printf("%f ", temp[j]);
136 }
137 printf("\n");
138 free(temp);
139 }
140 printf("\n");
141
142 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matA));
143 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matB));
144 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matC));
145 HIPSPARSE_CHECK(hipsparseDestroy(handle));
146
147 HIP_CHECK(hipFree(dA));
148 HIP_CHECK(hipFree(dB));
149 HIP_CHECK(hipFree(dcsr_row_ptrC));
150 HIP_CHECK(hipFree(dcsr_col_indC));
151 HIP_CHECK(hipFree(dcsr_valC));
152 HIP_CHECK(hipFree(dbuffer));
153
154 return 0;
155}
hipsparseSpSV_createDescr()#
-
hipsparseStatus_t hipsparseSpSV_createDescr(hipsparseSpSVDescr_t *descr)#
hipsparseSpSV_createDescrcreates a sparse matrix triangular solve descriptor. It should be destroyed at the end using hipsparseSpSV_destroyDescr().
hipsparseSpSV_destroyDescr()#
-
hipsparseStatus_t hipsparseSpSV_destroyDescr(hipsparseSpSVDescr_t descr)#
hipsparseSpSV_destroyDescrdestroys a sparse matrix triangular solve descriptor and releases all resources used by the descriptor.
hipsparseSpSV_bufferSize()#
-
hipsparseStatus_t hipsparseSpSV_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t opA, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstDnVecDescr_t x, const hipsparseDnVecDescr_t y, hipDataType computeType, hipsparseSpSVAlg_t alg, hipsparseSpSVDescr_t spsvDescr, size_t *pBufferSizeInBytes)#
hipsparseSpSV_bufferSizecomputes the size of the required user allocated buffer needed when solving the triangular linear system:\[ op(A) \cdot y := \alpha \cdot x, \]where \(op(A)\) is a sparse matrix in CSR or COO storage format, \(x\) and \(y\) are dense vectors.hipsparseSpSV_bufferSizesupports multiple combinations of data types and compute types. See hipsparseSpSV_solve for a complete listing of all the data type and compute type combinations available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] matrix operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] sparse matrix descriptor.
x – [in] dense vector descriptor.
y – [inout] dense vector descriptor.
computeType – [in] floating point precision for the SpSV computation.
alg – [in] SpSV algorithm for the SpSV computation.
spsvDescr – [in] SpSV descriptor.
pBufferSizeInBytes – [out] number of bytes of the temporary storage buffer.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_NOT_INITIALIZED –
handleis not initialized.HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,matA,x,y,spsvDescrorpBufferSizeInBytesis nullptr, oropAis invalid.HIPSPARSE_STATUS_NOT_SUPPORTED –
opA,computeTypeoralgis currently not supported.
hipsparseSpSV_analysis()#
-
hipsparseStatus_t hipsparseSpSV_analysis(hipsparseHandle_t handle, hipsparseOperation_t opA, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstDnVecDescr_t x, const hipsparseDnVecDescr_t y, hipDataType computeType, hipsparseSpSVAlg_t alg, hipsparseSpSVDescr_t spsvDescr, void *externalBuffer)#
hipsparseSpSV_analysisperforms the required analysis needed when solving the triangular linear system:\[ op(A) \cdot y := \alpha \cdot x, \]where \(op(A)\) is a sparse matrix in CSR or COO storage format, \(x\) and \(y\) are dense vectors.hipsparseSpSV_analysissupports multiple combinations of data types and compute types. See hipsparseSpSV_solve for a complete listing of all the data type and compute type combinations available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] matrix operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] matrix descriptor.
x – [in] vector descriptor.
y – [inout] vector descriptor.
computeType – [in] floating point precision for the SpSV computation.
alg – [in] SpSV algorithm for the SpSV computation.
spsvDescr – [in] SpSV descriptor.
externalBuffer – [out] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,matA,x,y,spsvDescrorexternalBufferpointer is invalid.HIPSPARSE_STATUS_NOT_SUPPORTED –
opA,computeTypeoralgis currently not supported.
hipsparseSpSV_solve()#
-
hipsparseStatus_t hipsparseSpSV_solve(hipsparseHandle_t handle, hipsparseOperation_t opA, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstDnVecDescr_t x, const hipsparseDnVecDescr_t y, hipDataType computeType, hipsparseSpSVAlg_t alg, hipsparseSpSVDescr_t spsvDescr)#
Sparse triangular solve.
hipsparseSpSV_solvesolves a triangular linear system of equations defined by a sparse \(m \times m\) square matrix \(op(A)\), given in CSR or COO storage format, 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}\]and where \(y\) is the dense solution vector and \(x\) is the dense right-hand side vector.Performing the above operation requires three steps. First, hipsparseSpSV_bufferSize must be called which will determine the size of the required temporary storage buffer. The user then allocates this buffer and calls hipsparseSpSV_analysis which will perform analysis on the sparse matrix \(op(A)\). Finally, the user completes the computation by calling
hipsparseSpSV_solve. The buffer size and preprecess routines only need to be called once for a given sparse matrix \(op(A)\) while the computation can be repeatedly used with different \(x\) and \(y\) vectors. Once all calls tohipsparseSpSV_solveare complete, the temporary buffer can be deallocated.hipsparseSpSV_solvesupports HIPSPARSE_INDEX_32I and HIPSPARSE_INDEX_64I index types for storing the row pointer and column indices arrays of the sparse matrices.hipsparseSpSV_solvesupports the following data types for \(op(A)\), \(x\), \(y\) and compute types for \(\alpha\):- Uniform Precisions:
A / X / Y / compute_type
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] matrix operation type.
alpha – [in] scalar \(\alpha\).
matA – [in] matrix descriptor.
x – [in] vector descriptor.
y – [inout] vector descriptor.
computeType – [in] floating point precision for the SpSV computation.
alg – [in] SpSV algorithm for the SpSV computation.
spsvDescr – [in] SpSV descriptor.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,matA,x,y, orspsvDescrpointer is invalid.HIPSPARSE_STATUS_NOT_SUPPORTED –
opA,computeTypeoralgis currently not supported.
1int main(int argc, char* argv[])
2{
3 // 1 0 0 0
4 // A = 4 2 0 0
5 // 0 3 7 0
6 // 0 0 0 1
7 int m = 4;
8
9 std::vector<int> hcsr_row_ptr = {0, 1, 3, 5, 6};
10 std::vector<int> hcsr_col_ind = {0, 0, 1, 1, 2, 3};
11 std::vector<float> hcsr_val = {1, 4, 2, 3, 7, 1};
12 std::vector<float> hx(m, 1.0f);
13 std::vector<float> hy(m, 0.0f);
14
15 // Scalar alpha
16 float alpha = 1.0f;
17
18 int nnz = hcsr_row_ptr[m] - hcsr_row_ptr[0];
19
20 // Offload data to device
21 int* dcsr_row_ptr;
22 int* dcsr_col_ind;
23 float* dcsr_val;
24 float* dx;
25 float* dy;
26 HIP_CHECK(hipMalloc((void**)&dcsr_row_ptr, sizeof(int) * (m + 1)));
27 HIP_CHECK(hipMalloc((void**)&dcsr_col_ind, sizeof(int) * nnz));
28 HIP_CHECK(hipMalloc((void**)&dcsr_val, sizeof(float) * nnz));
29 HIP_CHECK(hipMalloc((void**)&dx, sizeof(float) * m));
30 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * m));
31
32 HIP_CHECK(
33 hipMemcpy(dcsr_row_ptr, hcsr_row_ptr.data(), sizeof(int) * (m + 1), hipMemcpyHostToDevice));
34 HIP_CHECK(
35 hipMemcpy(dcsr_col_ind, hcsr_col_ind.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
36 HIP_CHECK(hipMemcpy(dcsr_val, hcsr_val.data(), sizeof(float) * nnz, hipMemcpyHostToDevice));
37 HIP_CHECK(hipMemcpy(dx, hx.data(), sizeof(float) * m, hipMemcpyHostToDevice));
38
39 hipsparseHandle_t handle;
40 hipsparseSpMatDescr_t matA;
41 hipsparseDnVecDescr_t vecX;
42 hipsparseDnVecDescr_t vecY;
43
44 hipsparseIndexType_t row_idx_type = HIPSPARSE_INDEX_32I;
45 hipsparseIndexType_t col_idx_type = HIPSPARSE_INDEX_32I;
46 hipDataType data_type = HIP_R_32F;
47 hipDataType computeType = HIP_R_32F;
48 hipsparseIndexBase_t idx_base = HIPSPARSE_INDEX_BASE_ZERO;
49 hipsparseOperation_t trans = HIPSPARSE_OPERATION_NON_TRANSPOSE;
50
51 HIPSPARSE_CHECK(hipsparseCreate(&handle));
52
53 // Create sparse matrix A
54 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
55 m,
56 m,
57 nnz,
58 dcsr_row_ptr,
59 dcsr_col_ind,
60 dcsr_val,
61 row_idx_type,
62 col_idx_type,
63 idx_base,
64 data_type));
65
66 // Create dense vector X
67 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecX, m, dx, data_type));
68
69 // Create dense vector Y
70 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, m, dy, data_type));
71
72 hipsparseSpSVDescr_t descr;
73 HIPSPARSE_CHECK(hipsparseSpSV_createDescr(&descr));
74
75 // Call spsv to get buffer size
76 size_t buffer_size;
77 HIPSPARSE_CHECK(hipsparseSpSV_bufferSize(handle,
78 trans,
79 &alpha,
80 matA,
81 vecX,
82 vecY,
83 computeType,
84 HIPSPARSE_SPSV_ALG_DEFAULT,
85 descr,
86 &buffer_size));
87
88 void* temp_buffer;
89 HIP_CHECK(hipMalloc((void**)&temp_buffer, buffer_size));
90
91 // Call spsv to perform analysis
92 HIPSPARSE_CHECK(hipsparseSpSV_analysis(handle,
93 trans,
94 &alpha,
95 matA,
96 vecX,
97 vecY,
98 computeType,
99 HIPSPARSE_SPSV_ALG_DEFAULT,
100 descr,
101 temp_buffer));
102
103 // Call spsv to perform computation
104 HIPSPARSE_CHECK(hipsparseSpSV_solve(
105 handle, trans, &alpha, matA, vecX, vecY, computeType, HIPSPARSE_SPSV_ALG_DEFAULT, descr));
106
107 // Copy result back to host
108 HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(float) * m, hipMemcpyDeviceToHost));
109
110 // Clear hipSPARSE
111 HIPSPARSE_CHECK(hipsparseSpSV_destroyDescr(descr));
112 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matA));
113 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecX));
114 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
115 HIPSPARSE_CHECK(hipsparseDestroy(handle));
116
117 // Clear device memory
118 HIP_CHECK(hipFree(dcsr_row_ptr));
119 HIP_CHECK(hipFree(dcsr_col_ind));
120 HIP_CHECK(hipFree(dcsr_val));
121 HIP_CHECK(hipFree(dx));
122 HIP_CHECK(hipFree(dy));
123 HIP_CHECK(hipFree(temp_buffer));
124
125 return 0;
126}
1int main(int argc, char* argv[])
2{
3 // 1 0 0 0
4 // A = 4 2 0 0
5 // 0 3 7 0
6 // 0 0 0 1
7 int m = 4;
8
9 int hcsr_row_ptr[] = {0, 1, 3, 5, 6};
10 int hcsr_col_ind[] = {0, 0, 1, 1, 2, 3};
11 float hcsr_val[] = {1, 4, 2, 3, 7, 1};
12 float* hx = (float*)malloc(m * sizeof(float));
13 float* hy = (float*)malloc(m * sizeof(float));
14
15 // Scalar alpha
16 float alpha = 1.0;
17
18 int nnz = hcsr_row_ptr[m] - hcsr_row_ptr[0];
19
20 // Offload data to device
21 int* dcsr_row_ptr;
22 int* dcsr_col_ind;
23 float* dcsr_val;
24 float* dx;
25 float* dy;
26 HIP_CHECK(hipMalloc((void**)&dcsr_row_ptr, sizeof(int) * (m + 1)));
27 HIP_CHECK(hipMalloc((void**)&dcsr_col_ind, sizeof(int) * nnz));
28 HIP_CHECK(hipMalloc((void**)&dcsr_val, sizeof(float) * nnz));
29 HIP_CHECK(hipMalloc((void**)&dx, sizeof(float) * m));
30 HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * m));
31
32 HIP_CHECK(hipMemcpy(dcsr_row_ptr, hcsr_row_ptr, sizeof(int) * (m + 1), hipMemcpyHostToDevice));
33 HIP_CHECK(hipMemcpy(dcsr_col_ind, hcsr_col_ind, sizeof(int) * nnz, hipMemcpyHostToDevice));
34 HIP_CHECK(hipMemcpy(dcsr_val, hcsr_val, sizeof(float) * nnz, hipMemcpyHostToDevice));
35 HIP_CHECK(hipMemcpy(dx, hx, sizeof(float) * m, hipMemcpyHostToDevice));
36
37 hipsparseHandle_t handle;
38 hipsparseSpMatDescr_t matA;
39 hipsparseDnVecDescr_t vecX;
40 hipsparseDnVecDescr_t vecY;
41
42 hipsparseIndexType_t row_idx_type = HIPSPARSE_INDEX_32I;
43 hipsparseIndexType_t col_idx_type = HIPSPARSE_INDEX_32I;
44 hipDataType data_type = HIP_R_32F;
45 hipDataType computeType = HIP_R_32F;
46 hipsparseIndexBase_t idx_base = HIPSPARSE_INDEX_BASE_ZERO;
47 hipsparseOperation_t trans = HIPSPARSE_OPERATION_NON_TRANSPOSE;
48
49 HIPSPARSE_CHECK(hipsparseCreate(&handle));
50
51 // Create sparse matrix A
52 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
53 m,
54 m,
55 nnz,
56 dcsr_row_ptr,
57 dcsr_col_ind,
58 dcsr_val,
59 row_idx_type,
60 col_idx_type,
61 idx_base,
62 data_type));
63
64 // Create dense vector X
65 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecX, m, dx, data_type));
66
67 // Create dense vector Y
68 HIPSPARSE_CHECK(hipsparseCreateDnVec(&vecY, m, dy, data_type));
69
70 hipsparseSpSVDescr_t descr;
71 HIPSPARSE_CHECK(hipsparseSpSV_createDescr(&descr));
72
73 // Call spsv to get buffer size
74 size_t buffer_size;
75 HIPSPARSE_CHECK(hipsparseSpSV_bufferSize(handle,
76 trans,
77 &alpha,
78 matA,
79 vecX,
80 vecY,
81 computeType,
82 HIPSPARSE_SPSV_ALG_DEFAULT,
83 descr,
84 &buffer_size));
85
86 void* temp_buffer;
87 HIP_CHECK(hipMalloc((void**)&temp_buffer, buffer_size));
88
89 // Call spsv to perform analysis
90 HIPSPARSE_CHECK(hipsparseSpSV_analysis(handle,
91 trans,
92 &alpha,
93 matA,
94 vecX,
95 vecY,
96 computeType,
97 HIPSPARSE_SPSV_ALG_DEFAULT,
98 descr,
99 temp_buffer));
100
101 // Call spsv to perform computation
102 HIPSPARSE_CHECK(hipsparseSpSV_solve(
103 handle, trans, &alpha, matA, vecX, vecY, computeType, HIPSPARSE_SPSV_ALG_DEFAULT, descr));
104
105 // Copy result back to host
106 HIP_CHECK(hipMemcpy(hy, dy, sizeof(float) * m, hipMemcpyDeviceToHost));
107
108 // Clear hipSPARSE
109 HIPSPARSE_CHECK(hipsparseSpSV_destroyDescr(descr));
110 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matA));
111 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecX));
112 HIPSPARSE_CHECK(hipsparseDestroyDnVec(vecY));
113 HIPSPARSE_CHECK(hipsparseDestroy(handle));
114
115 // Clear device memory
116 HIP_CHECK(hipFree(dcsr_row_ptr));
117 HIP_CHECK(hipFree(dcsr_col_ind));
118 HIP_CHECK(hipFree(dcsr_val));
119 HIP_CHECK(hipFree(dx));
120 HIP_CHECK(hipFree(dy));
121 HIP_CHECK(hipFree(temp_buffer));
122
123 return 0;
124}
hipsparseSpSM_createDescr()#
-
hipsparseStatus_t hipsparseSpSM_createDescr(hipsparseSpSMDescr_t *descr)#
hipsparseSpSM_createDescrcreates a sparse matrix triangular solve with multiple rhs descriptor. It should be destroyed at the end using hipsparseSpSM_destroyDescr().
hipsparseSpSM_destroyDescr()#
-
hipsparseStatus_t hipsparseSpSM_destroyDescr(hipsparseSpSMDescr_t descr)#
hipsparseSpSM_destroyDescrdestroys a sparse matrix triangular solve with multiple rhs descriptor and releases all resources used by the descriptor.
hipsparseSpSM_bufferSize()#
-
hipsparseStatus_t hipsparseSpSM_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstDnMatDescr_t matB, const hipsparseDnMatDescr_t matC, hipDataType computeType, hipsparseSpSMAlg_t alg, hipsparseSpSMDescr_t spsmDescr, size_t *pBufferSizeInBytes)#
hipsparseSpSM_bufferSizecomputes the size of the required user allocated buffer needed when solving the triangular linear system:\[ op(A) \cdot C := \alpha \cdot op(B), \]where \(op(A)\) is a square sparse matrix in CSR or COO storage format, \(B\) and \(C\) are dense matrices.hipsparseSpSM_bufferSizesupports multiple combinations of data types and compute types. See hipsparseSpSM_solve for a complete listing of all the data type and compute type combinations available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] matrix operation type for the sparse matrix \(A\).
opB – [in] matrix operation type for the dense matrix \(B\).
alpha – [in] scalar \(\alpha\).
matA – [in] sparse matrix descriptor.
matB – [in] dense matrix descriptor.
matC – [inout] dense matrix descriptor.
computeType – [in] floating point precision for the SpSM computation.
alg – [in] SpSM algorithm for the SpSM computation.
spsmDescr – [in] SpSM descriptor.
pBufferSizeInBytes – [out] number of bytes of the temporary storage buffer.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_NOT_INITIALIZED –
handleis not initialized.HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,matA,matB,matC,spsmDescrorpBufferSizeInBytesis nullptr, oropAoropBis invalid.HIPSPARSE_STATUS_NOT_SUPPORTED –
opA,opB,computeTypeoralgis currently not supported.
hipsparseSpSM_analysis()#
-
hipsparseStatus_t hipsparseSpSM_analysis(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstDnMatDescr_t matB, const hipsparseDnMatDescr_t matC, hipDataType computeType, hipsparseSpSMAlg_t alg, hipsparseSpSMDescr_t spsmDescr, void *externalBuffer)#
hipsparseSpSM_analysisperforms the required analysis needed when solving the triangular linear system:\[ op(A) \cdot C := \alpha \cdot op(B), \]where \(A\) is a sparse matrix in CSR or COO storage format, \(B\) and \(C\) are dense vectors.hipsparseSpSM_bufferSizesupports multiple combinations of data types and compute types. See hipsparseSpSM_solve for a complete listing of all the data type and compute type combinations available.- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] matrix operation type for the sparse matrix \(A\).
opB – [in] matrix operation type for the dense matrix \(B\).
alpha – [in] scalar \(\alpha\).
matA – [in] sparse matrix descriptor.
matB – [in] dense matrix descriptor.
matC – [inout] dense matrix descriptor.
computeType – [in] floating point precision for the SpSM computation.
alg – [in] SpSM algorithm for the SpSM computation.
spsmDescr – [in] SpSM descriptor.
externalBuffer – [out] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,matA,matB,matC,spsmDescrorexternalBufferpointer is invalid.HIPSPARSE_STATUS_NOT_SUPPORTED –
opA,opB,computeTypeoralgis currently not supported.
hipsparseSpSM_solve()#
-
hipsparseStatus_t hipsparseSpSM_solve(hipsparseHandle_t handle, hipsparseOperation_t opA, hipsparseOperation_t opB, const void *alpha, hipsparseConstSpMatDescr_t matA, hipsparseConstDnMatDescr_t matB, const hipsparseDnMatDescr_t matC, hipDataType computeType, hipsparseSpSMAlg_t alg, hipsparseSpSMDescr_t spsmDescr, void *externalBuffer)#
Sparse triangular system solve.
hipsparseSpSM_solvesolves a triangular linear system of equations defined by a sparse \(m \times m\) square matrix \(op(A)\), given in CSR or COO storage format, such that\[ op(A) \cdot C = \alpha \cdot op(B), \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if opA == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ A^T, & \text{if opB == HIPSPARSE_OPERATION_TRANSPOSE} \end{array} \right. \end{split}\]and\[\begin{split} op(B) = \left\{ \begin{array}{ll} B, & \text{if opA == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ B^T, & \text{if opB == HIPSPARSE_OPERATION_TRANSPOSE} \end{array} \right. \end{split}\]and where \(C\) is the dense solution matrix and \(B\) is the dense right-hand side matrix. Both \(B\) and \(C\) can be in row or column order.Performing the above operation requires three steps. First, the user calls hipsparseSpSM_bufferSize in order to determine the size of the required temporary storage buffer. The user then allocates this buffer and calls hipsparseSpSM_analysis which will perform analysis on the sparse matrix \(op(A)\). Finally, the user completes the computation by calling
hipsparseSpSM_solve. The buffer size and analysis routines only need to be called once for a given sparse matrix \(op(A)\) while the computation can be called repeatedly with different \(B\) and \(C\) matrices. Once all calls tohipsparseSpSM_solveare complete, the temporary buffer can be deallocated.As noted above, both \(B\) and \(C\) can be in row or column order (this includes mixing the order so that \(B\) is row order and \(C\) is column order and vice versa). When running on an AMD system with the rocSPARSE backend, the kernels solve the system assuming the matrices \(B\) and \(C\) are in row order as this provides the best memory access. This means that if the matrix \(C\) is not in row order and/or the matrix \(B\) is not row order (or \(B^{T}\) is not column order as this is equivalent to being in row order), then internally memory copies and/or transposing of data may be performed to get them into the correct order (possbily using extra buffer size). Once computation is completed, additional memory copies and/or transposing of data may be performed to get them back into the user arrays. For best performance and smallest required temporary storage buffers on an AMD system, use row order for the matrix \(C\) and row order for the matrix \(B\) (or column order if \(B\) is being transposed).
hipsparseSpSM_solvesupports HIPSPARSE_INDEX_32I and HIPSPARSE_INDEX_64I index precisions for storing the row pointer and column indices arrays of the sparse matrices.hipsparseSpSM_solvesupports the following data types for \(op(A)\), \(op(B)\), \(C\) and compute types for \(\alpha\):- Uniform Precisions:
A / B / C / compute_type
HIP_R_32F
HIP_R_64F
HIP_C_32F
HIP_C_64F
- Parameters:
handle – [in] handle to the hipsparse library context queue.
opA – [in] matrix operation type for the sparse matrix \(A\).
opB – [in] matrix operation type for the dense matrix \(B\).
alpha – [in] scalar \(\alpha\).
matA – [in] sparse matrix descriptor.
matB – [in] dense matrix descriptor.
matC – [inout] dense matrix descriptor.
computeType – [in] floating point precision for the SpSM computation.
alg – [in] SpSM algorithm for the SpSM computation.
spsmDescr – [in] SpSM descriptor.
externalBuffer – [out] temporary storage buffer allocated by the user.
- Return values:
HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.
HIPSPARSE_STATUS_INVALID_VALUE –
handle,alpha,matA,matB,matC,spsmDescrorexternalBufferpointer is invalid.HIPSPARSE_STATUS_NOT_SUPPORTED –
opA,opB,computeTypeoralgis currently not supported.
1int main(int argc, char* argv[])
2{
3 // 1 0 0 0
4 // A = 4 2 0 0
5 // 0 3 7 0
6 // 0 0 0 1
7 int m = 4;
8 int n = 2;
9
10 std::vector<int> hcsr_row_ptr = {0, 1, 3, 5, 6};
11 std::vector<int> hcsr_col_ind = {0, 0, 1, 1, 2, 3};
12 std::vector<float> hcsr_val = {1, 4, 2, 3, 7, 1};
13 std::vector<float> hB(m * n);
14 std::vector<float> hC(m * n);
15
16 for(int i = 0; i < n; i++)
17 {
18 for(int j = 0; j < m; j++)
19 {
20 hB[m * i + j] = static_cast<float>(i + 1);
21 }
22 }
23
24 // Scalar alpha
25 float alpha = 1.0f;
26
27 int nnz = hcsr_row_ptr[m] - hcsr_row_ptr[0];
28
29 // Offload data to device
30 int* dcsr_row_ptr;
31 int* dcsr_col_ind;
32 float* dcsr_val;
33 float* dB;
34 float* dC;
35 HIP_CHECK(hipMalloc((void**)&dcsr_row_ptr, sizeof(int) * (m + 1)));
36 HIP_CHECK(hipMalloc((void**)&dcsr_col_ind, sizeof(int) * nnz));
37 HIP_CHECK(hipMalloc((void**)&dcsr_val, sizeof(float) * nnz));
38 HIP_CHECK(hipMalloc((void**)&dB, sizeof(float) * m * n));
39 HIP_CHECK(hipMalloc((void**)&dC, sizeof(float) * m * n));
40
41 HIP_CHECK(
42 hipMemcpy(dcsr_row_ptr, hcsr_row_ptr.data(), sizeof(int) * (m + 1), hipMemcpyHostToDevice));
43 HIP_CHECK(
44 hipMemcpy(dcsr_col_ind, hcsr_col_ind.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
45 HIP_CHECK(hipMemcpy(dcsr_val, hcsr_val.data(), sizeof(float) * nnz, hipMemcpyHostToDevice));
46 HIP_CHECK(hipMemcpy(dB, hB.data(), sizeof(float) * m * n, hipMemcpyHostToDevice));
47
48 hipsparseHandle_t handle;
49 hipsparseSpMatDescr_t matA;
50 hipsparseDnMatDescr_t matB;
51 hipsparseDnMatDescr_t matC;
52
53 hipsparseIndexType_t row_idx_type = HIPSPARSE_INDEX_32I;
54 hipsparseIndexType_t col_idx_type = HIPSPARSE_INDEX_32I;
55 hipDataType dataType = HIP_R_32F;
56 hipDataType computeType = HIP_R_32F;
57 hipsparseIndexBase_t idxBase = HIPSPARSE_INDEX_BASE_ZERO;
58 hipsparseOperation_t transA = HIPSPARSE_OPERATION_NON_TRANSPOSE;
59 hipsparseOperation_t transB = HIPSPARSE_OPERATION_NON_TRANSPOSE;
60
61 HIPSPARSE_CHECK(hipsparseCreate(&handle));
62
63 // Create sparse matrix A
64 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
65 m,
66 m,
67 nnz,
68 dcsr_row_ptr,
69 dcsr_col_ind,
70 dcsr_val,
71 row_idx_type,
72 col_idx_type,
73 idxBase,
74 dataType));
75
76 // Create dense matrix B
77 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matB, m, n, m, dB, dataType, HIPSPARSE_ORDER_COL));
78
79 // Create dense matrix C
80 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matC, m, n, m, dC, dataType, HIPSPARSE_ORDER_COL));
81
82 hipsparseSpSMDescr_t descr;
83 HIPSPARSE_CHECK(hipsparseSpSM_createDescr(&descr));
84
85 // Call SpSM to get buffer size
86 size_t buffer_size;
87 HIPSPARSE_CHECK(hipsparseSpSM_bufferSize(handle,
88 transA,
89 transB,
90 &alpha,
91 matA,
92 matB,
93 matC,
94 computeType,
95 HIPSPARSE_SPSM_ALG_DEFAULT,
96 descr,
97 &buffer_size));
98
99 void* temp_buffer;
100 HIP_CHECK(hipMalloc((void**)&temp_buffer, buffer_size));
101
102 // Call SpSM to perform analysis
103 HIPSPARSE_CHECK(hipsparseSpSM_analysis(handle,
104 transA,
105 transB,
106 &alpha,
107 matA,
108 matB,
109 matC,
110 computeType,
111 HIPSPARSE_SPSM_ALG_DEFAULT,
112 descr,
113 temp_buffer));
114
115 // Call SpSM to perform computation
116 HIPSPARSE_CHECK(hipsparseSpSM_solve(handle,
117 transA,
118 transB,
119 &alpha,
120 matA,
121 matB,
122 matC,
123 computeType,
124 HIPSPARSE_SPSM_ALG_DEFAULT,
125 descr,
126 temp_buffer));
127
128 // Copy result back to host
129 HIP_CHECK(hipMemcpy(hC.data(), dC, sizeof(float) * m * n, hipMemcpyDeviceToHost));
130
131 // Clear hipSPARSE
132 HIPSPARSE_CHECK(hipsparseSpSM_destroyDescr(descr));
133 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matA));
134 HIPSPARSE_CHECK(hipsparseDestroyDnMat(matB));
135 HIPSPARSE_CHECK(hipsparseDestroyDnMat(matC));
136 HIPSPARSE_CHECK(hipsparseDestroy(handle));
137
138 // Clear device memory
139 HIP_CHECK(hipFree(dcsr_row_ptr));
140 HIP_CHECK(hipFree(dcsr_col_ind));
141 HIP_CHECK(hipFree(dcsr_val));
142 HIP_CHECK(hipFree(dB));
143 HIP_CHECK(hipFree(dC));
144 HIP_CHECK(hipFree(temp_buffer));
145
146 return 0;
147}
1int main(int argc, char* argv[])
2{
3 // 1 0 0 0
4 // A = 4 2 0 0
5 // 0 3 7 0
6 // 0 0 0 1
7 int m = 4;
8 int n = 2;
9
10 int hcsr_row_ptr[] = {0, 1, 3, 5, 6};
11 int hcsr_col_ind[] = {0, 0, 1, 1, 2, 3};
12 float hcsr_val[] = {1, 4, 2, 3, 7, 1};
13 float* hB = (float*)malloc((m * n) * sizeof(float));
14 float* hC = (float*)malloc((m * n) * sizeof(float));
15
16 for(int i = 0; i < n; i++)
17 {
18 for(int j = 0; j < m; j++)
19 {
20 hB[m * i + j] = (float)(i + 1);
21 }
22 }
23
24 // Scalar alpha
25 float alpha = 1.0;
26
27 int nnz = hcsr_row_ptr[m] - hcsr_row_ptr[0];
28
29 // Offload data to device
30 int* dcsr_row_ptr;
31 int* dcsr_col_ind;
32 float* dcsr_val;
33 float* dB;
34 float* dC;
35 HIP_CHECK(hipMalloc((void**)&dcsr_row_ptr, sizeof(int) * (m + 1)));
36 HIP_CHECK(hipMalloc((void**)&dcsr_col_ind, sizeof(int) * nnz));
37 HIP_CHECK(hipMalloc((void**)&dcsr_val, sizeof(float) * nnz));
38 HIP_CHECK(hipMalloc((void**)&dB, sizeof(float) * m * n));
39 HIP_CHECK(hipMalloc((void**)&dC, sizeof(float) * m * n));
40
41 HIP_CHECK(hipMemcpy(dcsr_row_ptr, hcsr_row_ptr, sizeof(int) * (m + 1), hipMemcpyHostToDevice));
42 HIP_CHECK(hipMemcpy(dcsr_col_ind, hcsr_col_ind, sizeof(int) * nnz, hipMemcpyHostToDevice));
43 HIP_CHECK(hipMemcpy(dcsr_val, hcsr_val, sizeof(float) * nnz, hipMemcpyHostToDevice));
44 HIP_CHECK(hipMemcpy(dB, hB, sizeof(float) * m * n, hipMemcpyHostToDevice));
45
46 hipsparseHandle_t handle;
47 hipsparseSpMatDescr_t matA;
48 hipsparseDnMatDescr_t matB;
49 hipsparseDnMatDescr_t matC;
50
51 hipsparseIndexType_t row_idx_type = HIPSPARSE_INDEX_32I;
52 hipsparseIndexType_t col_idx_type = HIPSPARSE_INDEX_32I;
53 hipDataType dataType = HIP_R_32F;
54 hipDataType computeType = HIP_R_32F;
55 hipsparseIndexBase_t idxBase = HIPSPARSE_INDEX_BASE_ZERO;
56 hipsparseOperation_t transA = HIPSPARSE_OPERATION_NON_TRANSPOSE;
57 hipsparseOperation_t transB = HIPSPARSE_OPERATION_NON_TRANSPOSE;
58
59 HIPSPARSE_CHECK(hipsparseCreate(&handle));
60
61 // Create sparse matrix A
62 HIPSPARSE_CHECK(hipsparseCreateCsr(&matA,
63 m,
64 m,
65 nnz,
66 dcsr_row_ptr,
67 dcsr_col_ind,
68 dcsr_val,
69 row_idx_type,
70 col_idx_type,
71 idxBase,
72 dataType));
73
74 // Create dense matrix B
75 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matB, m, n, m, dB, dataType, HIPSPARSE_ORDER_COL));
76
77 // Create dense matrix C
78 HIPSPARSE_CHECK(hipsparseCreateDnMat(&matC, m, n, m, dC, dataType, HIPSPARSE_ORDER_COL));
79
80 hipsparseSpSMDescr_t descr;
81 HIPSPARSE_CHECK(hipsparseSpSM_createDescr(&descr));
82
83 // Call SpSM to get buffer size
84 size_t buffer_size;
85 HIPSPARSE_CHECK(hipsparseSpSM_bufferSize(handle,
86 transA,
87 transB,
88 &alpha,
89 matA,
90 matB,
91 matC,
92 computeType,
93 HIPSPARSE_SPSM_ALG_DEFAULT,
94 descr,
95 &buffer_size));
96
97 void* temp_buffer;
98 HIP_CHECK(hipMalloc((void**)&temp_buffer, buffer_size));
99
100 // Call SpSM to perform analysis
101 HIPSPARSE_CHECK(hipsparseSpSM_analysis(handle,
102 transA,
103 transB,
104 &alpha,
105 matA,
106 matB,
107 matC,
108 computeType,
109 HIPSPARSE_SPSM_ALG_DEFAULT,
110 descr,
111 temp_buffer));
112
113 // Call SpSM to perform computation
114 HIPSPARSE_CHECK(hipsparseSpSM_solve(handle,
115 transA,
116 transB,
117 &alpha,
118 matA,
119 matB,
120 matC,
121 computeType,
122 HIPSPARSE_SPSM_ALG_DEFAULT,
123 descr,
124 temp_buffer));
125
126 // Copy result back to host
127 HIP_CHECK(hipMemcpy(hC, dC, sizeof(float) * m * n, hipMemcpyDeviceToHost));
128
129 // Clear hipSPARSE
130 HIPSPARSE_CHECK(hipsparseSpSM_destroyDescr(descr));
131 HIPSPARSE_CHECK(hipsparseDestroyMatDescr(matA));
132 HIPSPARSE_CHECK(hipsparseDestroyDnMat(matB));
133 HIPSPARSE_CHECK(hipsparseDestroyDnMat(matC));
134 HIPSPARSE_CHECK(hipsparseDestroy(handle));
135
136 // Clear device memory
137 HIP_CHECK(hipFree(dcsr_row_ptr));
138 HIP_CHECK(hipFree(dcsr_col_ind));
139 HIP_CHECK(hipFree(dcsr_val));
140 HIP_CHECK(hipFree(dB));
141 HIP_CHECK(hipFree(dC));
142 HIP_CHECK(hipFree(temp_buffer));
143
144 return 0;
145}