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[x_ind[i]] = y[x_ind[i]] + alpha * x_val[i]; }
- Example
// Number of non-zeros of the sparse vector int nnz = 3; // Sparse index vector int hx_ind[3] = {0, 3, 5}; // Sparse value vector double hx_val[3] = {1.0, 2.0, 3.0}; // Dense vector double hy[9] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0}; // Scalar alpha double alpha = 3.7; // Index base hipsparseIndexBase_t idx_base = HIPSPARSE_INDEX_BASE_ZERO; // Offload data to device int* dx_ind; double* dx_val; double* dy; hipMalloc((void**)&dx_ind, sizeof(int) * nnz); hipMalloc((void**)&dx_val, sizeof(double) * nnz); hipMalloc((void**)&dy, sizeof(double) * 9); hipMemcpy(dx_ind, hx_ind, sizeof(int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dx_val, hx_val, sizeof(double) * nnz, hipMemcpyHostToDevice); hipMemcpy(dy, hy, sizeof(double) * 9, hipMemcpyHostToDevice); // hipSPARSE handle hipsparseHandle_t handle; hipsparseCreate(&handle); // Call daxpyi to perform y = y + alpha * x hipsparseDaxpyi(handle, nnz, &alpha, dx_val, dx_ind, dy, idx_base); // Copy result back to host hipMemcpy(hy, dy, sizeof(double) * 9, hipMemcpyDeviceToHost); // Clear hipSPARSE hipsparseDestroy(handle); // Clear device memory hipFree(dx_ind); hipFree(dx_val); hipFree(dy);
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
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\[ \text{result} := y^T x \]for(i = 0; i < nnz; ++i) { result += x_val[i] * y[x_ind[i]]; }
- Example
// Number of non-zeros of the sparse vector int nnz = 3; // Sparse index vector int hx_ind[3] = {0, 3, 5}; // Sparse value vector float hx_val[3] = {1.0f, 2.0f, 3.0f}; // Dense vector float hy[9] = {1.0f, 2.0f, 3.0f, 4.0f, 5.0f, 6.0f, 7.0f, 8.0f, 9.0f}; // Index base hipsparseIndexBase_t idx_base = HIPSPARSE_INDEX_BASE_ZERO; // Offload data to device int* dx_ind; float* dx_val; float* dy; hipMalloc((void**)&dx_ind, sizeof(int) * nnz); hipMalloc((void**)&dx_val, sizeof(float) * nnz); hipMalloc((void**)&dy, sizeof(float) * 9); hipMemcpy(dx_ind, hx_ind, sizeof(int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dx_val, hx_val, sizeof(float) * nnz, hipMemcpyHostToDevice); hipMemcpy(dy, hy, sizeof(float) * 9, hipMemcpyHostToDevice); // hipSPARSE handle hipsparseHandle_t handle; hipsparseCreate(&handle); // Call sdoti to compute the dot product float dot; hipsparseSdoti(handle, nnz, dx_val, dx_ind, dy, &dot, idx_base); // Clear hipSPARSE hipsparseDestroy(handle); // Clear device memory hipFree(dx_ind); hipFree(dx_val); hipFree(dy);
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
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\[ \text{result} := \bar{x}^H y \]for(i = 0; i < nnz; ++i) { result += conj(x_val[i]) * y[x_ind[i]]; }
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
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 inx_ind
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]]; }
- Example
// Number of non-zeros of the sparse vector int nnz = 3; // Sparse index vector int hx_ind[3] = {0, 3, 5}; // Sparse value vector float hx_val[3]; // Dense vector float hy[9] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0}; // Index base hipsparseIndexBase_t idx_base = HIPSPARSE_INDEX_BASE_ZERO; // Offload data to device int* dx_ind; float* dx_val; float* dy; hipMalloc((void**)&dx_ind, sizeof(int) * nnz); hipMalloc((void**)&dx_val, sizeof(float) * nnz); hipMalloc((void**)&dy, sizeof(float) * 9); hipMemcpy(dx_ind, hx_ind, sizeof(int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dy, hy, sizeof(float) * 9, hipMemcpyHostToDevice); // hipSPARSE handle hipsparseHandle_t handle; hipsparseCreate(&handle); // Call sgthr hipsparseSgthr(handle, nnz, dy, dx_val, dx_ind, idx_base); // Copy result back to host hipMemcpy(hx_val, dx_val, sizeof(float) * nnz, hipMemcpyDeviceToHost); // Clear hipSPARSE hipsparseDestroy(handle); // Clear device memory hipFree(dx_ind); hipFree(dx_val); hipFree(dy);
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
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 inx_ind
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) { x_val[i] = y[x_ind[i]]; y[x_ind[i]] = 0; }
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
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 = 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; }
- Example
// Number of non-zeros of the sparse vector int nnz = 3; // Sparse index vector int hx_ind[3] = {0, 3, 5}; // Sparse value vector float hx_val[3] = {1.0f, 2.0f, 3.0f}; // Dense vector float hy[9] = {1.0f, 2.0f, 3.0f, 4.0f, 5.0f, 6.0f, 7.0f, 8.0f, 9.0f}; // c and s float c = 3.7; float s = 1.3; // Index base hipsparseIndexBase_t idx_base = HIPSPARSE_INDEX_BASE_ZERO; // Offload data to device int* dx_ind; float* dx_val; float* dy; hipMalloc((void**)&dx_ind, sizeof(int) * nnz); hipMalloc((void**)&dx_val, sizeof(float) * nnz); hipMalloc((void**)&dy, sizeof(float) * 9); hipMemcpy(dx_ind, hx_ind, sizeof(int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dx_val, hx_val, sizeof(float) * nnz, hipMemcpyHostToDevice); hipMemcpy(dy, hy, sizeof(float) * 9, hipMemcpyHostToDevice); // hipSPARSE handle hipsparseHandle_t handle; hipsparseCreate(&handle); // Call sroti hipsparseSroti(handle, nnz, dx_val, dx_ind, dy, &c, &s, idx_base); // Copy result back to host hipMemcpy(hx_val, dx_val, sizeof(float) * nnz, hipMemcpyDeviceToHost); hipMemcpy(hy, dy, sizeof(float) * 9, hipMemcpyDeviceToHost); // Clear hipSPARSE hipsparseDestroy(handle); // Clear device memory hipFree(dx_ind); hipFree(dx_val); hipFree(dy);
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
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 inx_ind
from the sparse vector \(x\) into the dense vector \(y\). Indices of \(y\) that are not listed inx_ind
remain unchanged.for(i = 0; i < nnz; ++i) { y[x_ind[i]] = x_val[i]; }
- Example
// Number of non-zeros of the sparse vector int nnz = 3; // Sparse index vector int hx_ind[3] = {0, 3, 5}; // Sparse value vector float hx_val[3] = {9.0, 2.0, 3.0}; // Dense vector float hy[9] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0}; // Index base hipsparseIndexBase_t idx_base = HIPSPARSE_INDEX_BASE_ZERO; // Offload data to device int* dx_ind; float* dx_val; float* dy; hipMalloc((void**)&dx_ind, sizeof(int) * nnz); hipMalloc((void**)&dx_val, sizeof(float) * nnz); hipMalloc((void**)&dy, sizeof(float) * 9); hipMemcpy(dx_ind, hx_ind, sizeof(int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dx_val, hx_val, sizeof(float) * nnz, hipMemcpyHostToDevice); hipMemcpy(dy, hy, sizeof(float) * 9, hipMemcpyHostToDevice); // hipSPARSE handle hipsparseHandle_t handle; hipsparseCreate(&handle); // Call ssctr hipsparseSsctr(handle, nnz, dx_val, dx_ind, dy, idx_base); // Copy result back to host hipMemcpy(hy, dy, sizeof(float) * 9, hipMemcpyDeviceToHost); // Clear hipSPARSE hipsparseDestroy(handle); // Clear device memory hipFree(dx_ind); hipFree(dx_val); hipFree(dy);
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.