Sparse level 1 functions#

The sparse level 1 routines describe operations between a vector in sparse format and a vector in dense format. This section describes all hipSPARSE level 1 sparse linear algebra functions.

hipsparseXaxpyi()#

hipsparseStatus_t hipsparseSaxpyi(hipsparseHandle_t handle, int nnz, const float *alpha, const float *xVal, const int *xInd, float *y, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseDaxpyi(hipsparseHandle_t handle, int nnz, const double *alpha, const double *xVal, const int *xInd, double *y, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseCaxpyi(hipsparseHandle_t handle, int nnz, const hipComplex *alpha, const hipComplex *xVal, const int *xInd, hipComplex *y, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseZaxpyi(hipsparseHandle_t handle, int nnz, const hipDoubleComplex *alpha, const hipDoubleComplex *xVal, const int *xInd, hipDoubleComplex *y, hipsparseIndexBase_t idxBase)#

Scale a sparse vector and add it to a dense vector.

hipsparseXaxpyi multiplies the sparse vector \(x\) with scalar \(\alpha\) and adds the result to the dense vector \(y\), such that

\[ y := y + \alpha \cdot x \]

for(i = 0; i < nnz; ++i)
{
    y[xInd[i]] = y[xInd[i]] + alpha * xVal[i];
}

Deprecated:

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

Note

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

Note

If nnz is zero, the function returns successfully without modifying y. Duplicate indices in xInd will result in the corresponding values being added multiple times to the same location in y.

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

  • nnz[in] number of non-zero entries of vector \(x\). Must be non-negative.

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

  • xVal[in] array of nnz elements containing the values of \(x\).

  • xInd[in] array of nnz elements containing the indices of the non-zero values of \(x\).

  • y[inout] array of values in dense format. Must be pre-allocated with sufficient size to accommodate all indices specified in xInd.

  • idxBase[in] index base. HIPSPARSE_INDEX_BASE_ZERO for zero-based indexing or HIPSPARSE_INDEX_BASE_ONE for one-based indexing.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle is nullptr, nnz is negative, alpha, xVal, xInd or y is nullptr when nnz is greater than zero, or idxBase is neither HIPSPARSE_INDEX_BASE_ZERO nor HIPSPARSE_INDEX_BASE_ONE.

 1int main(int argc, char* argv[])
 2{
 3    // Number of non-zeros of the sparse vector
 4    const int nnz = 3;
 5
 6    // Number of entries in the dense vector
 7    const int size = 9;
 8
 9    // Sparse index vector
10    std::vector<int> hxInd = {0, 3, 5};
11
12    // Sparse value vector
13    std::vector<double> hxVal = {1.0, 2.0, 3.0};
14
15    // Dense vector
16    std::vector<double> hy = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0};
17
18    // Scalar alpha
19    double alpha = 3.7;
20
21    // Index base
22    hipsparseIndexBase_t idxBase = HIPSPARSE_INDEX_BASE_ZERO;
23
24    // Offload data to device
25    int*    dxInd;
26    double* dxVal;
27    double* dy;
28
29    HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
30    HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(double) * nnz));
31    HIP_CHECK(hipMalloc((void**)&dy, sizeof(double) * size));
32
33    HIP_CHECK(hipMemcpy(dxInd, hxInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
34    HIP_CHECK(hipMemcpy(dxVal, hxVal.data(), sizeof(double) * nnz, hipMemcpyHostToDevice));
35    HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(double) * size, hipMemcpyHostToDevice));
36
37    // hipSPARSE handle
38    hipsparseHandle_t handle;
39    HIPSPARSE_CHECK(hipsparseCreate(&handle));
40
41    // Call daxpyi to perform y = y + alpha * x
42    HIPSPARSE_CHECK(hipsparseDaxpyi(handle, nnz, &alpha, dxVal, dxInd, dy, idxBase));
43
44    // Copy result back to host
45    HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(double) * size, hipMemcpyDeviceToHost));
46
47    std::cout << "hy" << std::endl;
48    for(int i = 0; i < size; i++)
49    {
50        std::cout << hy[i] << " ";
51    }
52    std::cout << std::endl;
53
54    // Clear hipSPARSE
55    HIPSPARSE_CHECK(hipsparseDestroy(handle));
56
57    // Clear device memory
58    HIP_CHECK(hipFree(dxInd));
59    HIP_CHECK(hipFree(dxVal));
60    HIP_CHECK(hipFree(dy));
61    return 0;
62}

hipsparseXdoti()#

hipsparseStatus_t hipsparseSdoti(hipsparseHandle_t handle, int nnz, const float *xVal, const int *xInd, const float *y, float *result, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseDdoti(hipsparseHandle_t handle, int nnz, const double *xVal, const int *xInd, const double *y, double *result, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseCdoti(hipsparseHandle_t handle, int nnz, const hipComplex *xVal, const int *xInd, const hipComplex *y, hipComplex *result, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseZdoti(hipsparseHandle_t handle, int nnz, const hipDoubleComplex *xVal, const int *xInd, const hipDoubleComplex *y, hipDoubleComplex *result, hipsparseIndexBase_t idxBase)#

Compute the dot product of a sparse vector with a dense vector.

hipsparseXdoti computes the dot product of the sparse vector \(x\) with the dense vector \(y\), such that

\[ result := y^T x \]

result = 0
for(i = 0; i < nnz; ++i)
{
    result += xVal[i] * y[xInd[i]];
}

Deprecated:

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

Note

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

Note

If nnz is zero, the function returns successfully with result set to zero.

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

  • nnz[in] number of non-zero entries of vector \(x\). Must be non-negative.

  • xVal[in] array of nnz values containing the elements of \(x\).

  • xInd[in] array of nnz elements containing the indices of the non-zero values of \(x\).

  • y[in] array of values in dense format. Must be pre-allocated with sufficient size to accommodate all indices specified in xInd.

  • result[out] pointer to the result, can be host or device memory.

  • idxBase[in] index base. HIPSPARSE_INDEX_BASE_ZERO for zero-based indexing or HIPSPARSE_INDEX_BASE_ONE for one-based indexing.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle or result is nullptr, nnz is negative, xVal, xInd or y is nullptr when nnz is greater than zero, or idxBase is neither HIPSPARSE_INDEX_BASE_ZERO nor HIPSPARSE_INDEX_BASE_ONE.

  • HIPSPARSE_STATUS_ALLOC_FAILED – the buffer for the dot product reduction could not be allocated.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

 1int main(int argc, char* argv[])
 2{
 3    // Number of non-zeros of the sparse vector
 4    const int nnz = 3;
 5
 6    // Number of entries in the dense vector
 7    const 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    // Index base
19    hipsparseIndexBase_t idxBase = HIPSPARSE_INDEX_BASE_ZERO;
20
21    // Offload data to device
22    int*   dxInd;
23    float* dxVal;
24    float* dy;
25
26    HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
27    HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
28    HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
29
30    HIP_CHECK(hipMemcpy(dxInd, hxInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
31    HIP_CHECK(hipMemcpy(dxVal, hxVal.data(), sizeof(float) * nnz, hipMemcpyHostToDevice));
32    HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(float) * size, hipMemcpyHostToDevice));
33
34    // hipSPARSE handle
35    hipsparseHandle_t handle;
36    HIPSPARSE_CHECK(hipsparseCreate(&handle));
37
38    // Call sdoti to compute the dot product
39    float dot;
40    HIPSPARSE_CHECK(hipsparseSdoti(handle, nnz, dxVal, dxInd, dy, &dot, idxBase));
41
42    std::cout << "dot: " << dot << std::endl;
43
44    // Clear hipSPARSE
45    HIPSPARSE_CHECK(hipsparseDestroy(handle));
46
47    // Clear device memory
48    HIP_CHECK(hipFree(dxInd));
49    HIP_CHECK(hipFree(dxVal));
50    HIP_CHECK(hipFree(dy));
51
52    return 0;
53}

hipsparseXdotci()#

hipsparseStatus_t hipsparseCdotci(hipsparseHandle_t handle, int nnz, const hipComplex *xVal, const int *xInd, const hipComplex *y, hipComplex *result, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseZdotci(hipsparseHandle_t handle, int nnz, const hipDoubleComplex *xVal, const int *xInd, const hipDoubleComplex *y, hipDoubleComplex *result, hipsparseIndexBase_t idxBase)#

Compute the dot product of a complex conjugate sparse vector with a dense vector.

hipsparseXdotci computes the dot product of the complex conjugate sparse vector \(x\) with the dense vector \(y\), such that

\[ result := \bar{x}^H y \]

result = 0
for(i = 0; i < nnz; ++i)
{
    result += conj(xVal[i]) * y[xInd[i]];
}

Deprecated:

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

Note

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

Note

If nnz is zero, the function returns successfully with result set to zero.

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

  • nnz[in] number of non-zero entries of vector \(x\). Must be non-negative.

  • xVal[in] array of nnz values containing the elements of \(x\).

  • xInd[in] array of nnz elements containing the indices of the non-zero values of \(x\).

  • y[in] array of values in dense format. Must be pre-allocated with sufficient size to accommodate all indices specified in xInd.

  • result[out] pointer to the result, can be host or device memory.

  • idxBase[in] index base. HIPSPARSE_INDEX_BASE_ZERO for zero-based indexing or HIPSPARSE_INDEX_BASE_ONE for one-based indexing.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle or result is nullptr, nnz is negative, xVal, xInd or y is nullptr when nnz is greater than zero, or idxBase is neither HIPSPARSE_INDEX_BASE_ZERO nor HIPSPARSE_INDEX_BASE_ONE.

  • HIPSPARSE_STATUS_ALLOC_FAILED – the buffer for the dot product reduction could not be allocated.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXgthr()#

hipsparseStatus_t hipsparseSgthr(hipsparseHandle_t handle, int nnz, const float *y, float *xVal, const int *xInd, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseDgthr(hipsparseHandle_t handle, int nnz, const double *y, double *xVal, const int *xInd, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseCgthr(hipsparseHandle_t handle, int nnz, const hipComplex *y, hipComplex *xVal, const int *xInd, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseZgthr(hipsparseHandle_t handle, int nnz, const hipDoubleComplex *y, hipDoubleComplex *xVal, const int *xInd, hipsparseIndexBase_t idxBase)#

Gather elements from a dense vector and store them into a sparse vector.

hipsparseXgthr gathers the elements that are listed in xInd from the dense vector \(y\) and stores them in the sparse vector \(x\).

for(i = 0; i < nnz; ++i)
{
    xVal[i] = y[xInd[i]];
}

Deprecated:

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

Note

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

Note

If nnz is zero, the function returns successfully without modifying xVal.

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

  • nnz[in] number of non-zero entries of \(x\). Must be non-negative.

  • y[in] array of values in dense format. Must be pre-allocated with sufficient size to accommodate all indices specified in xInd.

  • xVal[out] array of nnz elements that will contain the gathered values of \(x\).

  • xInd[in] array of nnz elements containing the indices of the non-zero values of \(x\).

  • idxBase[in] index base. HIPSPARSE_INDEX_BASE_ZERO for zero-based indexing or HIPSPARSE_INDEX_BASE_ONE for one-based indexing.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle is nullptr, nnz is negative, y, xVal or xInd is nullptr when nnz is greater than zero, or idxBase is neither HIPSPARSE_INDEX_BASE_ZERO nor HIPSPARSE_INDEX_BASE_ONE.

 1int main(int argc, char* argv[])
 2{
 3    // Number of non-zeros of the sparse vector
 4    const int nnz = 3;
 5
 6    // Number of entries in the dense vector
 7    const 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(nnz);
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    // Index base
19    hipsparseIndexBase_t idxBase = HIPSPARSE_INDEX_BASE_ZERO;
20
21    // Offload data to device
22    int*   dxInd;
23    float* dxVal;
24    float* dy;
25
26    HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
27    HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
28    HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
29
30    HIP_CHECK(hipMemcpy(dxInd, hxInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
31    HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(float) * size, hipMemcpyHostToDevice));
32
33    // hipSPARSE handle
34    hipsparseHandle_t handle;
35    HIPSPARSE_CHECK(hipsparseCreate(&handle));
36
37    // Call sgthr
38    HIPSPARSE_CHECK(hipsparseSgthr(handle, nnz, dy, dxVal, dxInd, idxBase));
39
40    // Copy result back to host
41    HIP_CHECK(hipMemcpy(hxVal.data(), dxVal, sizeof(float) * nnz, hipMemcpyDeviceToHost));
42
43    std::cout << "hxVal" << std::endl;
44    for(int i = 0; i < nnz; i++)
45    {
46        std::cout << hxVal[i] << " ";
47    }
48    std::cout << std::endl;
49
50    // Clear hipSPARSE
51    HIPSPARSE_CHECK(hipsparseDestroy(handle));
52
53    // Clear device memory
54    HIP_CHECK(hipFree(dxInd));
55    HIP_CHECK(hipFree(dxVal));
56    HIP_CHECK(hipFree(dy));
57
58    return 0;
59}

hipsparseXgthrz()#

hipsparseStatus_t hipsparseSgthrz(hipsparseHandle_t handle, int nnz, float *y, float *xVal, const int *xInd, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseDgthrz(hipsparseHandle_t handle, int nnz, double *y, double *xVal, const int *xInd, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseCgthrz(hipsparseHandle_t handle, int nnz, hipComplex *y, hipComplex *xVal, const int *xInd, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseZgthrz(hipsparseHandle_t handle, int nnz, hipDoubleComplex *y, hipDoubleComplex *xVal, const int *xInd, hipsparseIndexBase_t idxBase)#

Gather and zero out elements from a dense vector and store them into a sparse vector.

hipsparseXgthrz gathers the elements that are listed in xInd from the dense vector \(y\) and stores them in the sparse vector \(x\). The gathered elements in \(y\) are replaced by zero.

for(i = 0; i < nnz; ++i)
{
    xVal[i]    = y[xInd[i]];
    y[xInd[i]] = 0;
}

Deprecated:

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

Note

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

Note

If nnz is zero, the function returns successfully without modifying xVal or y.

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

  • nnz[in] number of non-zero entries of \(x\). Must be non-negative.

  • y[inout] array of values in dense format. Must be pre-allocated with sufficient size to accommodate all indices specified in xInd. Gathered elements are set to zero.

  • xVal[out] array of nnz elements that will contain the gathered values of \(x\).

  • xInd[in] array of nnz elements containing the indices of the non-zero values of \(x\).

  • idxBase[in] index base. HIPSPARSE_INDEX_BASE_ZERO for zero-based indexing or HIPSPARSE_INDEX_BASE_ONE for one-based indexing.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle is nullptr, nnz is negative, y, xVal or xInd is nullptr when nnz is greater than zero, or idxBase is neither HIPSPARSE_INDEX_BASE_ZERO nor HIPSPARSE_INDEX_BASE_ONE.

hipsparseXroti()#

hipsparseStatus_t hipsparseSroti(hipsparseHandle_t handle, int nnz, float *xVal, const int *xInd, float *y, const float *c, const float *s, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseDroti(hipsparseHandle_t handle, int nnz, double *xVal, const int *xInd, double *y, const double *c, const double *s, hipsparseIndexBase_t idxBase)#

Apply Givens rotation to a dense and a sparse vector.

hipsparseXroti applies 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 = xVal[i];
    y_tmp = y[xInd[i]];

    xVal[i]    = c * x_tmp + s * y_tmp;
    y[xInd[i]] = c * y_tmp - s * x_tmp;
}

Deprecated:

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

Note

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

Note

If nnz is zero, the function returns successfully without modifying xVal or y.

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

  • nnz[in] number of non-zero entries of \(x\). Must be non-negative.

  • xVal[inout] array of nnz elements containing the non-zero values of \(x\).

  • xInd[in] array of nnz elements containing the indices of the non-zero values of \(x\).

  • y[inout] array of values in dense format. Must be pre-allocated with sufficient size to accommodate all indices specified in xInd.

  • c[in] pointer to the cosine element of \(G\), can be on host or device.

  • s[in] pointer to the sine element of \(G\), can be on host or device.

  • idxBase[in] index base. HIPSPARSE_INDEX_BASE_ZERO for zero-based indexing or HIPSPARSE_INDEX_BASE_ONE for one-based indexing.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, c or s is nullptr, nnz is negative, xVal, xInd or y is nullptr when nnz is greater than zero, or idxBase is neither HIPSPARSE_INDEX_BASE_ZERO nor HIPSPARSE_INDEX_BASE_ONE.

 1int main(int argc, char* argv[])
 2{
 3    // Number of non-zeros of the sparse vector
 4    const int nnz = 3;
 5
 6    // Number of entries in the dense vector
 7    const 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    // c and s
19    float c = 3.7;
20    float s = 1.3;
21
22    // Index base
23    hipsparseIndexBase_t idxBase = HIPSPARSE_INDEX_BASE_ZERO;
24
25    // Offload data to device
26    int*   dxInd;
27    float* dxVal;
28    float* dy;
29
30    HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
31    HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
32    HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
33
34    HIP_CHECK(hipMemcpy(dxInd, hxInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
35    HIP_CHECK(hipMemcpy(dxVal, hxVal.data(), sizeof(float) * nnz, hipMemcpyHostToDevice));
36    HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(float) * size, hipMemcpyHostToDevice));
37
38    // hipSPARSE handle
39    hipsparseHandle_t handle;
40    HIPSPARSE_CHECK(hipsparseCreate(&handle));
41
42    // Call sroti
43    HIPSPARSE_CHECK(hipsparseSroti(handle, nnz, dxVal, dxInd, dy, &c, &s, idxBase));
44
45    // Copy result back to host
46    HIP_CHECK(hipMemcpy(hxVal.data(), dxVal, sizeof(float) * nnz, hipMemcpyDeviceToHost));
47    HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(float) * size, hipMemcpyDeviceToHost));
48
49    std::cout << "hxVal" << std::endl;
50    for(int i = 0; i < nnz; i++)
51    {
52        std::cout << hxVal[i] << " ";
53    }
54    std::cout << std::endl;
55
56    std::cout << "hy" << std::endl;
57    for(int i = 0; i < size; i++)
58    {
59        std::cout << hy[i] << " ";
60    }
61    std::cout << std::endl;
62
63    // Clear hipSPARSE
64    HIPSPARSE_CHECK(hipsparseDestroy(handle));
65
66    // Clear device memory
67    HIP_CHECK(hipFree(dxInd));
68    HIP_CHECK(hipFree(dxVal));
69    HIP_CHECK(hipFree(dy));
70
71    return 0;
72}

hipsparseXsctr()#

hipsparseStatus_t hipsparseSsctr(hipsparseHandle_t handle, int nnz, const float *xVal, const int *xInd, float *y, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseDsctr(hipsparseHandle_t handle, int nnz, const double *xVal, const int *xInd, double *y, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseCsctr(hipsparseHandle_t handle, int nnz, const hipComplex *xVal, const int *xInd, hipComplex *y, hipsparseIndexBase_t idxBase)#
hipsparseStatus_t hipsparseZsctr(hipsparseHandle_t handle, int nnz, const hipDoubleComplex *xVal, const int *xInd, hipDoubleComplex *y, hipsparseIndexBase_t idxBase)#

Scatter elements from a dense vector across a sparse vector.

hipsparseXsctr scatters the elements that are listed in xInd from the sparse vector \(x\) into the dense vector \(y\). Indices of \(y\) that are not listed in xInd remain unchanged.

for(i = 0; i < nnz; ++i)
{
    y[xInd[i]] = xVal[i];
}

Deprecated:

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

Note

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

Note

If nnz is zero, the function returns successfully without modifying y. Duplicate indices in xInd will result in the last value being written to y.

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

  • nnz[in] number of non-zero entries of \(x\). Must be non-negative.

  • xVal[in] array of nnz elements containing the non-zero values of \(x\).

  • xInd[in] array of nnz elements containing the indices of the non-zero values of \(x\).

  • y[inout] array of values in dense format. Must be pre-allocated with sufficient size to accommodate all indices specified in xInd.

  • idxBase[in] index base. HIPSPARSE_INDEX_BASE_ZERO for zero-based indexing or HIPSPARSE_INDEX_BASE_ONE for one-based indexing.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_NOT_INITIALIZEDhandle is not initialized.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle is nullptr, nnz is negative, xVal, xInd or y is nullptr when nnz is greater than zero, or idxBase is neither HIPSPARSE_INDEX_BASE_ZERO nor HIPSPARSE_INDEX_BASE_ONE.

 1int main(int argc, char* argv[])
 2{
 3    // Number of non-zeros of the sparse vector
 4    const int nnz = 3;
 5
 6    // Number of entries in the dense vector
 7    const 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 = {9.0, 2.0, 3.0};
14
15    // Dense vector
16    std::vector<float> hy = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0};
17
18    // Index base
19    hipsparseIndexBase_t idxBase = HIPSPARSE_INDEX_BASE_ZERO;
20
21    // Offload data to device
22    int*   dxInd;
23    float* dxVal;
24    float* dy;
25
26    HIP_CHECK(hipMalloc((void**)&dxInd, sizeof(int) * nnz));
27    HIP_CHECK(hipMalloc((void**)&dxVal, sizeof(float) * nnz));
28    HIP_CHECK(hipMalloc((void**)&dy, sizeof(float) * size));
29
30    HIP_CHECK(hipMemcpy(dxInd, hxInd.data(), sizeof(int) * nnz, hipMemcpyHostToDevice));
31    HIP_CHECK(hipMemcpy(dxVal, hxVal.data(), sizeof(float) * nnz, hipMemcpyHostToDevice));
32    HIP_CHECK(hipMemcpy(dy, hy.data(), sizeof(float) * size, hipMemcpyHostToDevice));
33
34    // hipSPARSE handle
35    hipsparseHandle_t handle;
36    HIPSPARSE_CHECK(hipsparseCreate(&handle));
37
38    // Call ssctr
39    HIPSPARSE_CHECK(hipsparseSsctr(handle, nnz, dxVal, dxInd, dy, idxBase));
40
41    // Copy result back to host
42    HIP_CHECK(hipMemcpy(hy.data(), dy, sizeof(float) * size, hipMemcpyDeviceToHost));
43
44    std::cout << "hy" << std::endl;
45    for(int i = 0; i < size; i++)
46    {
47        std::cout << hy[i] << " ";
48    }
49    std::cout << std::endl;
50
51    // Clear hipSPARSE
52    HIPSPARSE_CHECK(hipsparseDestroy(handle));
53
54    // Clear device memory
55    HIP_CHECK(hipFree(dxInd));
56    HIP_CHECK(hipFree(dxVal));
57    HIP_CHECK(hipFree(dy));
58
59    return 0;
60}