Sparse Level 2 Functions#
This module holds all sparse level 2 routines.
The sparse level 2 routines describe operations between a matrix in sparse format and a vector in dense format.
hipsparseXcsrmv()#
-
hipsparseStatus_t hipsparseScsrmv(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, const float *alpha, const hipsparseMatDescr_t descrA, const float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const float *x, const float *beta, float *y)#
-
hipsparseStatus_t hipsparseDcsrmv(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, const double *alpha, const hipsparseMatDescr_t descrA, const double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const double *x, const double *beta, double *y)#
-
hipsparseStatus_t hipsparseCcsrmv(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, const hipComplex *alpha, const hipsparseMatDescr_t descrA, const hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const hipComplex *x, const hipComplex *beta, hipComplex *y)#
-
hipsparseStatus_t hipsparseZcsrmv(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, const hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, const hipDoubleComplex *x, const hipDoubleComplex *beta, hipDoubleComplex *y)#
Sparse matrix vector multiplication using CSR storage format.
hipsparseXcsrmv
multiplies the scalar \(\alpha\) with a sparse \(m \times n\) matrix, defined in CSR storage format, and the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if 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}\]for(i = 0; i < m; ++i) { y[i] = beta * y[i]; for(j = csrRowPtr[i]; j < csrRowPtr[i + 1]; ++j) { y[i] = y[i] + alpha * csrVal[j] * x[csrColInd[j]]; } }
- Example
// hipSPARSE handle hipsparseHandle_t handle; hipsparseCreate(&handle); // alpha * ( 1.0 0.0 2.0 ) * ( 1.0 ) + beta * ( 4.0 ) = ( 31.1 ) // ( 3.0 0.0 4.0 ) * ( 2.0 ) ( 5.0 ) = ( 62.0 ) // ( 5.0 6.0 0.0 ) * ( 3.0 ) ( 6.0 ) = ( 70.7 ) // ( 7.0 0.0 8.0 ) * ( 7.0 ) = ( 123.8 ) int m = 4; int n = 3; int nnz = 8; // CSR row pointers int hcsrRowPtr[5] = {0, 2, 4, 6, 8}; // CSR column indices int hcsrColInd[8] = {0, 2, 0, 2, 0, 1, 0, 2}; // CSR values double hcsrVal[8] = {1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0}; // Transposition of the matrix hipsparseOperation_t trans = HIPSPARSE_OPERATION_NON_TRANSPOSE; // Scalar alpha and beta double alpha = 3.7; double beta = 1.3; // x and y double hx[3] = {1.0, 2.0, 3.0}; double hy[4] = {4.0, 5.0, 6.0, 7.0}; // Matrix descriptor hipsparseMatDescr_t descr; hipsparseCreateMatDescr(&descr); // Offload data to device int* dcsrRowPtr; int* dcsrColInd; double* dcsrVal; double* dx; double* dy; hipMalloc((void**)&dcsrRowPtr, sizeof(int) * (m + 1)); hipMalloc((void**)&dcsrColInd, sizeof(int) * nnz); hipMalloc((void**)&dcsrVal, sizeof(double) * nnz); hipMalloc((void**)&dx, sizeof(double) * n); hipMalloc((void**)&dy, sizeof(double) * m); hipMemcpy(dcsrRowPtr, hcsrRowPtr, sizeof(int) * (m + 1), hipMemcpyHostToDevice); hipMemcpy(dcsrColInd, hcsrColInd, sizeof(int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dcsrVal, hcsrVal, sizeof(double) * nnz, hipMemcpyHostToDevice); hipMemcpy(dx, hx, sizeof(double) * n, hipMemcpyHostToDevice); hipMemcpy(dy, hy, sizeof(double) * m, hipMemcpyHostToDevice); // Call dcsrmv to perform y = alpha * A x + beta * y hipsparseDcsrmv(handle, trans, m, n, nnz, &alpha, descr, dcsrVal, dcsrRowPtr, dcsrColInd, dx, &beta, dy); // Copy result back to host hipMemcpy(hy, dy, sizeof(double) * m, hipMemcpyDeviceToHost); // Clear hipSPARSE hipsparseDestroyMatDescr(descr); hipsparseDestroy(handle); // Clear device memory hipFree(dcsrRowPtr); hipFree(dcsrColInd); hipFree(dcsrVal); hipFree(dx); 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.
Note
Currently, only
trans
== HIPSPARSE_OPERATION_NON_TRANSPOSE is supported.
hipsparseXcsrsv2_zeroPivot()#
-
hipsparseStatus_t hipsparseXcsrsv2_zeroPivot(hipsparseHandle_t handle, csrsv2Info_t info, int *position)#
Sparse triangular solve using CSR storage format.
hipsparseXcsrsv2_zeroPivot
returns HIPSPARSE_STATUS_ZERO_PIVOT, if either a structural or numerical zero has been found during hipsparseScsrsv2_solve(), hipsparseDcsrsv2_solve(), hipsparseCcsrsv2_solve() or hipsparseZcsrsv2_solve() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored inposition
, using same index base as the CSR matrix.position
can be in host or device memory. If no zero pivot has been found,position
is set to -1 and HIPSPARSE_STATUS_SUCCESS is returned instead.Note
hipsparseXcsrsv2_zeroPivot
is a blocking function. It might influence performance negatively.
hipsparseXcsrsv2_bufferSize()#
-
hipsparseStatus_t hipsparseScsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, int *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseDcsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, int *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseCcsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, int *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseZcsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, int *pBufferSizeInBytes)#
Sparse triangular solve using CSR storage format.
hipsparseXcsrsv2_bufferSize
returns the size of the temporary storage buffer in bytes that is required by hipsparseScsrsv2_analysis(), hipsparseDcsrsv2_analysis(), hipsparseCcsrsv2_analysis(), hipsparseZcsrsv2_analysis(), hipsparseScsrsv2_solve(), hipsparseDcsrsv2_solve(), hipsparseCcsrsv2_solve() and hipsparseZcsrsv2_solve(). The temporary storage buffer must be allocated by the user.
hipsparseXcsrsv2_bufferSizeExt()#
-
hipsparseStatus_t hipsparseScsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, size_t *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseDcsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, size_t *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseCcsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, size_t *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseZcsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, size_t *pBufferSizeInBytes)#
Sparse triangular solve using CSR storage format.
hipsparseXcsrsv2_bufferSizeExt
returns the size of the temporary storage buffer in bytes that is required by hipsparseScsrsv2_analysis(), hipsparseDcsrsv2_analysis(), hipsparseCcsrsv2_analysis(), hipsparseZcsrsv2_analysis(), hipsparseScsrsv2_solve(), hipsparseDcsrsv2_solve(), hipsparseCcsrsv2_solve() and hipsparseZcsrsv2_solve(). The temporary storage buffer must be allocated by the user.
hipsparseXcsrsv2_analysis()#
-
hipsparseStatus_t hipsparseScsrsv2_analysis(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, const float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
-
hipsparseStatus_t hipsparseDcsrsv2_analysis(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, const double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
-
hipsparseStatus_t hipsparseCcsrsv2_analysis(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, const hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
-
hipsparseStatus_t hipsparseZcsrsv2_analysis(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipsparseMatDescr_t descrA, const hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
Sparse triangular solve using CSR storage format.
hipsparseXcsrsv2_analysis
performs the analysis step for hipsparseScsrsv2_solve(), hipsparseDcsrsv2_solve(), hipsparseCcsrsv2_solve() and hipsparseZcsrsv2_solve().Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
hipsparseXcsrsv2_solve()#
-
hipsparseStatus_t hipsparseScsrsv2_solve(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const float *alpha, const hipsparseMatDescr_t descrA, const float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, const float *f, float *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
-
hipsparseStatus_t hipsparseDcsrsv2_solve(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const double *alpha, const hipsparseMatDescr_t descrA, const double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, const double *f, double *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
-
hipsparseStatus_t hipsparseCcsrsv2_solve(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipComplex *alpha, const hipsparseMatDescr_t descrA, const hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, const hipComplex *f, hipComplex *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
-
hipsparseStatus_t hipsparseZcsrsv2_solve(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int nnz, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, const hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrsv2Info_t info, const hipDoubleComplex *f, hipDoubleComplex *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
Sparse triangular solve using CSR storage format.
hipsparseXcsrsv2_solve
solves a sparse triangular linear system of a sparse \(m \times m\) matrix, defined in CSR storage format, a dense solution vector \(y\) and the right-hand side \(x\) that is multiplied by \(\alpha\), such that\[ op(A) \cdot y = \alpha \cdot x, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ A^T, & \text{if trans == HIPSPARSE_OPERATION_TRANSPOSE} \\ A^H, & \text{if trans == HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE} \end{array} \right. \end{split}\]hipsparseXcsrsv2_solve
requires a user allocated temporary buffer. Its size is returned by hipsparseXcsrsv2_bufferSize() or hipsparseXcsrsv2_bufferSizeExt(). Furthermore, analysis meta data is required. It can be obtained by hipsparseXcsrsv2_analysis().hipsparseXcsrsv2_solve
reports the first zero pivot (either numerical or structural zero). The zero pivot status can be checked calling hipsparseXcsrsv2_zeroPivot(). If hipsparseDiagType_t == HIPSPARSE_DIAG_TYPE_UNIT, no zero pivot will be reported, even if \(A_{j,j} = 0\) for some \(j\).- Example
// hipSPARSE handle hipsparseHandle_t handle; hipsparseCreate(&handle); // alpha * ( 1.0 0.0 2.0 0.0 ) * ( x_0 ) = ( 32.0 ) // ( 3.0 2.0 4.0 1.0 ) * ( x_1 ) = ( 14.7 ) // ( 5.0 6.0 1.0 3.0 ) * ( x_2 ) = ( 33.6 ) // ( 7.0 0.0 8.0 0.6 ) * ( x_3 ) = ( 10.0 ) int m = 4; int nnz = 13; // CSR row pointers int hcsrRowPtr[5] = {0, 2, 6, 10, 13}; // CSR column indices int hcsrColInd[13] = {0, 2, 0, 1, 2, 3, 0, 1, 2, 3, 0, 2, 3}; // CSR values double hcsrVal[13] = {1.0, 2.0, 3.0, 2.0, 4.0, 1.0, 5.0, 6.0, 1.0, 3.0, 7.0, 8.0, 0.6}; // Transposition of the matrix hipsparseOperation_t trans = HIPSPARSE_OPERATION_NON_TRANSPOSE; hipsparseSolvePolicy_t policy = HIPSPARSE_SOLVE_POLICY_USE_LEVEL; // Scalar alpha double alpha = 1.0; // f and x double hf[4] = {32.0, 14.7, 33.6, 10.0}; double hx[4]; // Matrix descriptor hipsparseMatDescr_t descr; hipsparseCreateMatDescr(&descr); // Set index base on descriptor hipsparseSetMatIndexBase(descr, HIPSPARSE_INDEX_BASE_ZERO); // Set fill mode on descriptor hipsparseSetMatFillMode(descr, HIPSPARSE_FILL_MODE_LOWER); // Set diag type on descriptor hipsparseSetMatDiagType(descr, HIPSPARSE_DIAG_TYPE_UNIT); // Csrsv info csrsv2Info_t info; hipsparseCreateCsrsv2Info(&info); // Offload data to device int* dcsrRowPtr; int* dcsrColInd; double* dcsrVal; double* df; double* dx; hipMalloc((void**)&dcsrRowPtr, sizeof(int) * (m + 1)); hipMalloc((void**)&dcsrColInd, sizeof(int) * nnz); hipMalloc((void**)&dcsrVal, sizeof(double) * nnz); hipMalloc((void**)&df, sizeof(double) * m); hipMalloc((void**)&dx, sizeof(double) * m); hipMemcpy(dcsrRowPtr, hcsrRowPtr, sizeof(int) * (m + 1), hipMemcpyHostToDevice); hipMemcpy(dcsrColInd, hcsrColInd, sizeof(int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dcsrVal, hcsrVal, sizeof(double) * nnz, hipMemcpyHostToDevice); hipMemcpy(df, hf, sizeof(double) * m, hipMemcpyHostToDevice); int bufferSize = 0; hipsparseDcsrsv2_bufferSize(handle, trans, m, nnz, descr, dcsrVal, dcsrRowPtr, dcsrColInd, info, &bufferSize); void* dbuffer = nullptr; hipMalloc((void**)&dbuffer, bufferSize); hipsparseDcsrsv2_analysis(handle, trans, m, nnz, descr, dcsrVal, dcsrRowPtr, dcsrColInd, info, policy, dbuffer); // Call dcsrsv to perform alpha * A * x = f hipsparseDcsrsv2_solve(handle, trans, m, nnz, &alpha, descr, dcsrVal, dcsrRowPtr, dcsrColInd, info, df, dx, policy, dbuffer); // Copy result back to host hipMemcpy(hx, dx, sizeof(double) * m, hipMemcpyDeviceToHost); // Clear hipSPARSE hipsparseDestroyMatDescr(descr); hipsparseDestroyCsrsv2Info(info); hipsparseDestroy(handle); // Clear device memory hipFree(dcsrRowPtr); hipFree(dcsrColInd); hipFree(dcsrVal); hipFree(df); hipFree(dx); hipFree(dbuffer);
Note
The sparse CSR matrix has to be sorted. This can be achieved by calling hipsparseXcsrsort().
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
Currently, only
trans
== HIPSPARSE_OPERATION_NON_TRANSPOSE andtrans
== HIPSPARSE_OPERATION_TRANSPOSE is supported.
hipsparseXhybmv()#
-
hipsparseStatus_t hipsparseShybmv(hipsparseHandle_t handle, hipsparseOperation_t transA, const float *alpha, const hipsparseMatDescr_t descrA, const hipsparseHybMat_t hybA, const float *x, const float *beta, float *y)#
-
hipsparseStatus_t hipsparseDhybmv(hipsparseHandle_t handle, hipsparseOperation_t transA, const double *alpha, const hipsparseMatDescr_t descrA, const hipsparseHybMat_t hybA, const double *x, const double *beta, double *y)#
-
hipsparseStatus_t hipsparseChybmv(hipsparseHandle_t handle, hipsparseOperation_t transA, const hipComplex *alpha, const hipsparseMatDescr_t descrA, const hipsparseHybMat_t hybA, const hipComplex *x, const hipComplex *beta, hipComplex *y)#
-
hipsparseStatus_t hipsparseZhybmv(hipsparseHandle_t handle, hipsparseOperation_t transA, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, const hipsparseHybMat_t hybA, const hipDoubleComplex *x, const hipDoubleComplex *beta, hipDoubleComplex *y)#
Sparse matrix vector multiplication using HYB storage format.
hipsparseXhybmv
multiplies the scalar \(\alpha\) with a sparse \(m \times n\) matrix, defined in HYB storage format, and the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\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}\]- Example
// hipSPARSE handle hipsparseHandle_t handle; hipsparseCreate(&handle); // A sparse matrix // 1 0 3 4 // 0 0 5 1 // 0 2 0 0 // 4 0 0 8 int hAptr[5] = {0, 3, 5, 6, 8}; int hAcol[8] = {0, 2, 3, 2, 3, 1, 0, 3}; double hAval[8] = {1.0, 3.0, 4.0, 5.0, 1.0, 2.0, 4.0, 8.0}; int m = 4; int n = 4; int nnz = 8; double halpha = 1.0; double hbeta = 0.0; double hx[4] = {1.0, 2.0, 3.0, 4.0}; double hy[4] = {4.0, 5.0, 6.0, 7.0}; // Matrix descriptor hipsparseMatDescr_t descrA; hipsparseCreateMatDescr(&descrA); // Offload data to device int* dAptr = NULL; int* dAcol = NULL; double* dAval = NULL; double* dx = NULL; double* dy = NULL; hipMalloc((void**)&dAptr, sizeof(int) * (m + 1)); hipMalloc((void**)&dAcol, sizeof(int) * nnz); hipMalloc((void**)&dAval, sizeof(double) * nnz); hipMalloc((void**)&dx, sizeof(double) * n); hipMalloc((void**)&dy, sizeof(double) * m); hipMemcpy(dAptr, hAptr, sizeof(int) * (m + 1), hipMemcpyHostToDevice); hipMemcpy(dAcol, hAcol, sizeof(int) * nnz, hipMemcpyHostToDevice); hipMemcpy(dAval, hAval, sizeof(double) * nnz, hipMemcpyHostToDevice); hipMemcpy(dx, hx, sizeof(double) * n, hipMemcpyHostToDevice); // Convert CSR matrix to HYB format hipsparseHybMat_t hybA; hipsparseCreateHybMat(&hybA); hipsparseDcsr2hyb(handle, m, n, descrA, dAval, dAptr, dAcol, hybA, 0, HIPSPARSE_HYB_PARTITION_AUTO); // Clean up CSR structures hipFree(dAptr); hipFree(dAcol); hipFree(dAval); // Call hipsparse hybmv hipsparseDhybmv(handle, HIPSPARSE_OPERATION_NON_TRANSPOSE, &halpha, descrA, hybA, dx, &hbeta, dy); // Copy result back to host hipMemcpy(hy, dy, sizeof(double) * m, hipMemcpyDeviceToHost); // Clear up on device hipsparseDestroyHybMat(hybA); hipsparseDestroyMatDescr(descrA); hipsparseDestroy(handle); hipFree(dx); 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.
Note
Currently, only
trans
== HIPSPARSE_OPERATION_NON_TRANSPOSE is supported.
hipsparseXbsrmv()#
-
hipsparseStatus_t hipsparseSbsrmv(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nb, int nnzb, const float *alpha, const hipsparseMatDescr_t descrA, const float *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, const float *x, const float *beta, float *y)#
-
hipsparseStatus_t hipsparseDbsrmv(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nb, int nnzb, const double *alpha, const hipsparseMatDescr_t descrA, const double *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, const double *x, const double *beta, double *y)#
-
hipsparseStatus_t hipsparseCbsrmv(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nb, int nnzb, const hipComplex *alpha, const hipsparseMatDescr_t descrA, const hipComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, const hipComplex *x, const hipComplex *beta, hipComplex *y)#
-
hipsparseStatus_t hipsparseZbsrmv(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nb, int nnzb, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, const hipDoubleComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, const hipDoubleComplex *x, const hipDoubleComplex *beta, hipDoubleComplex *y)#
Sparse matrix vector multiplication using BSR storage format.
hipsparseXbsrmv
multiplies the scalar \(\alpha\) with a sparse \((mb \cdot \text{blockDim}) \times (nb \cdot \text{blockDim})\) matrix, defined in BSR storage format, and the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\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}\]- Example
// hipSPARSE handle hipsparseHandle_t handle; hipsparseCreate(&handle); // alpha * ( 1.0 0.0 2.0 ) * ( 1.0 ) + beta * ( 4.0 ) = ( 31.1 ) // ( 3.0 0.0 4.0 ) * ( 2.0 ) ( 5.0 ) = ( 62.0 ) // ( 5.0 6.0 0.0 ) * ( 3.0 ) ( 6.0 ) = ( 70.7 ) // ( 7.0 0.0 8.0 ) * ( 7.0 ) = ( 123.8 ) // BSR block dimension int bsr_dim = 2; // Number of block rows and columns int mb = 2; int nb = 2; // Number of non-zero blocks int nnzb = 4; // BSR row pointers int hbsrRowPtr[3] = {0, 2, 4}; // BSR column indices int hbsrColInd[4] = {0, 1, 0, 1}; // BSR values double hbsrVal[16] = {1.0, 3.0, 0.0, 0.0, 2.0, 4.0, 0.0, 0.0, 5.0, 7.0, 6.0, 0.0, 0.0, 8.0, 0.0, 0.0}; // Block storage in column major hipsparseDirection_t dir = HIPSPARSE_DIRECTION_COLUMN; // Transposition of the matrix hipsparseOperation_t trans = HIPSPARSE_OPERATION_NON_TRANSPOSE; // Scalar alpha and beta double alpha = 3.7; double beta = 1.3; // x and y double hx[4] = {1.0, 2.0, 3.0, 0.0}; double hy[4] = {4.0, 5.0, 6.0, 7.0}; // Matrix descriptor hipsparseMatDescr_t descr; hipsparseCreateMatDescr(&descr); // Offload data to device int* dbsrRowPtr; int* dbsrColInd; double* dbsrVal; double* dx; double* dy; hipMalloc((void**)&dbsrRowPtr, sizeof(int) * (mb + 1)); hipMalloc((void**)&dbsrColInd, sizeof(int) * nnzb); hipMalloc((void**)&dbsrVal, sizeof(double) * nnzb * bsr_dim * bsr_dim); hipMalloc((void**)&dx, sizeof(double) * nb * bsr_dim); hipMalloc((void**)&dy, sizeof(double) * mb * bsr_dim); hipMemcpy(dbsrRowPtr, hbsrRowPtr, sizeof(int) * (mb + 1), hipMemcpyHostToDevice); hipMemcpy(dbsrColInd, hbsrColInd, sizeof(int) * nnzb, hipMemcpyHostToDevice); hipMemcpy(dbsrVal, hbsrVal, sizeof(double) * nnzb * bsr_dim * bsr_dim, hipMemcpyHostToDevice); hipMemcpy(dx, hx, sizeof(double) * nb * bsr_dim, hipMemcpyHostToDevice); hipMemcpy(dy, hy, sizeof(double) * mb * bsr_dim, hipMemcpyHostToDevice); // Call dbsrmv to perform y = alpha * A x + beta * y hipsparseDbsrmv(handle, dir, trans, mb, nb, nnzb, &alpha, descr, dbsrVal, dbsrRowPtr, dbsrColInd, bsr_dim, dx, &beta, dy); // Copy result back to host hipMemcpy(hy, dy, sizeof(double) * mb * bsr_dim, hipMemcpyDeviceToHost); // Clear hipSPARSE hipsparseDestroyMatDescr(descr); hipsparseDestroy(handle); // Clear device memory hipFree(dbsrRowPtr); hipFree(dbsrColInd); hipFree(dbsrVal); hipFree(dx); 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.
Note
Currently, only
trans
== HIPSPARSE_OPERATION_NON_TRANSPOSE is supported.
hipsparseXbsrxmv()#
-
hipsparseStatus_t hipsparseSbsrxmv(hipsparseHandle_t handle, hipsparseDirection_t dir, hipsparseOperation_t trans, int sizeOfMask, int mb, int nb, int nnzb, const float *alpha, const hipsparseMatDescr_t descr, const float *bsrVal, const int *bsrMaskPtr, const int *bsrRowPtr, const int *bsrEndPtr, const int *bsrColInd, int blockDim, const float *x, const float *beta, float *y)#
-
hipsparseStatus_t hipsparseDbsrxmv(hipsparseHandle_t handle, hipsparseDirection_t dir, hipsparseOperation_t trans, int sizeOfMask, int mb, int nb, int nnzb, const double *alpha, const hipsparseMatDescr_t descr, const double *bsrVal, const int *bsrMaskPtr, const int *bsrRowPtr, const int *bsrEndPtr, const int *bsrColInd, int blockDim, const double *x, const double *beta, double *y)#
-
hipsparseStatus_t hipsparseCbsrxmv(hipsparseHandle_t handle, hipsparseDirection_t dir, hipsparseOperation_t trans, int sizeOfMask, int mb, int nb, int nnzb, const hipComplex *alpha, const hipsparseMatDescr_t descr, const hipComplex *bsrVal, const int *bsrMaskPtr, const int *bsrRowPtr, const int *bsrEndPtr, const int *bsrColInd, int blockDim, const hipComplex *x, const hipComplex *beta, hipComplex *y)#
-
hipsparseStatus_t hipsparseZbsrxmv(hipsparseHandle_t handle, hipsparseDirection_t dir, hipsparseOperation_t trans, int sizeOfMask, int mb, int nb, int nnzb, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descr, const hipDoubleComplex *bsrVal, const int *bsrMaskPtr, const int *bsrRowPtr, const int *bsrEndPtr, const int *bsrColInd, int blockDim, const hipDoubleComplex *x, const hipDoubleComplex *beta, hipDoubleComplex *y)#
Sparse matrix vector multiplication with mask operation using BSR storage format.
hipsparseXbsrxmv
multiplies the scalar \(\alpha\) with a sparse \((mb \cdot \text{blockDim}) \times (nb \cdot \text{blockDim})\) modified matrix, defined in BSR storage format, and the dense vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \left( \alpha \cdot op(A) \cdot x + \beta \cdot y \right)\left( \text{mask} \right), \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ A^T, & \text{if trans == HIPSPARSE_OPERATION_TRANSPOSE} \\ A^H, & \text{if trans == HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE} \end{array} \right. \end{split}\]The \(\text{mask}\) is defined as an array of block row indices. The input sparse matrix is defined with a modified BSR storage format where the beginning and the end of each row is defined with two arrays,
bsrRowPtr
andbsr_end_ptr
(both of sizemb
), rather the usualbsrRowPtr
of sizemb
+ 1.Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
Currently, only
trans
== HIPSPARSE_OPERATION_NON_TRANSPOSE is supported. Currently,blockDim
== 1 is not supported.
hipsparseXbsrsv2_zeroPivot()#
-
hipsparseStatus_t hipsparseXbsrsv2_zeroPivot(hipsparseHandle_t handle, bsrsv2Info_t info, int *position)#
Sparse triangular solve using BSR storage format.
hipsparseXbsrsv2_zeroPivot
returns HIPSPARSE_STATUS_ZERO_PIVOT, if either a structural or numerical zero has been found during hipsparseXbsrsv2_analysis() or hipsparseXbsrsv2_solve() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored inposition
, using same index base as the BSR matrix.position
can be in host or device memory. If no zero pivot has been found,position
is set to -1 and HIPSPARSE_STATUS_SUCCESS is returned instead.Note
hipsparseXbsrsv2_zeroPivot
is a blocking function. It might influence performance negatively.
hipsparseXbsrsv2_bufferSize()#
-
hipsparseStatus_t hipsparseSbsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, float *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, int *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseDbsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, double *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, int *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseCbsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, int *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseZbsrsv2_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipDoubleComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, int *pBufferSizeInBytes)#
Sparse triangular solve using BSR storage format.
hipsparseXbsrsv2_bufferSize
returns the size of the temporary storage buffer in bytes that is required by hipsparseXbsrsv2_analysis() and hipsparseXbsrsv2_solve(). The temporary storage buffer must be allocated by the user.
hipsparseXbsrsv2_bufferSizeExt()#
-
hipsparseStatus_t hipsparseSbsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, float *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, size_t *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseDbsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, double *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, size_t *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseCbsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, size_t *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseZbsrsv2_bufferSizeExt(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipDoubleComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, size_t *pBufferSizeInBytes)#
Sparse triangular solve using BSR storage format.
hipsparseXbsrsv2_bufferSizeExt
returns the size of the temporary storage buffer in bytes that is required by hipsparseXbsrsv2_analysis() and hipsparseXbsrsv2_solve(). The temporary storage buffer must be allocated by the user.
hipsparseXbsrsv2_analysis()#
-
hipsparseStatus_t hipsparseSbsrsv2_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, const float *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
-
hipsparseStatus_t hipsparseDbsrsv2_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, const double *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
-
hipsparseStatus_t hipsparseCbsrsv2_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, const hipComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
-
hipsparseStatus_t hipsparseZbsrsv2_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipsparseMatDescr_t descrA, const hipDoubleComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
Sparse triangular solve using BSR storage format.
hipsparseXbsrsv2_analysis
performs the analysis step for hipsparseXbsrsv2_solve().Note
If the matrix sparsity pattern changes, the gathered information will become invalid.
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
hipsparseXbsrsv2_solve()#
-
hipsparseStatus_t hipsparseSbsrsv2_solve(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const float *alpha, const hipsparseMatDescr_t descrA, const float *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, const float *f, float *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
-
hipsparseStatus_t hipsparseDbsrsv2_solve(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const double *alpha, const hipsparseMatDescr_t descrA, const double *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, const double *f, double *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
-
hipsparseStatus_t hipsparseCbsrsv2_solve(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipComplex *alpha, const hipsparseMatDescr_t descrA, const hipComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, const hipComplex *f, hipComplex *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
-
hipsparseStatus_t hipsparseZbsrsv2_solve(hipsparseHandle_t handle, hipsparseDirection_t dirA, hipsparseOperation_t transA, int mb, int nnzb, const hipDoubleComplex *alpha, const hipsparseMatDescr_t descrA, const hipDoubleComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrsv2Info_t info, const hipDoubleComplex *f, hipDoubleComplex *x, hipsparseSolvePolicy_t policy, void *pBuffer)#
Sparse triangular solve using BSR storage format.
hipsparseXbsrsv2_solve
solves a sparse triangular linear system of a sparse \(m \times m\) matrix, defined in BSR storage format, a dense solution vector \(y\) and the right-hand side \(x\) that is multiplied by \(\alpha\), such that\[ op(A) \cdot y = \alpha \cdot x, \]with\[\begin{split} op(A) = \left\{ \begin{array}{ll} A, & \text{if trans == HIPSPARSE_OPERATION_NON_TRANSPOSE} \\ A^T, & \text{if trans == HIPSPARSE_OPERATION_TRANSPOSE} \\ A^H, & \text{if trans == HIPSPARSE_OPERATION_CONJUGATE_TRANSPOSE} \end{array} \right. \end{split}\]hipsparseXbsrsv2_solve
requires a user allocated temporary buffer. Its size is returned by hipsparseXbsrsv2_bufferSize() or hipsparseXbsrsv2_bufferSizeExt(). Furthermore, analysis meta data is required. It can be obtained by hipsparseXbsrsv2_analysis().hipsparseXbsrsv2_solve
reports the first zero pivot (either numerical or structural zero). The zero pivot status can be checked calling hipsparseXbsrsv2_zeroPivot(). If hipsparseDiagType_t == HIPSPARSE_DIAG_TYPE_UNIT, no zero pivot will be reported, even if \(A_{j,j} = 0\) for some \(j\).- Example
// hipSPARSE handle hipsparseHandle_t handle; hipsparseCreate(&handle); // A = ( 1.0 0.0 0.0 0.0 ) // ( 2.0 3.0 0.0 0.0 ) // ( 4.0 5.0 6.0 0.0 ) // ( 7.0 0.0 8.0 9.0 ) // // with bsr_dim = 2 // // ------------------- // = | 1.0 0.0 | 0.0 0.0 | // | 2.0 3.0 | 0.0 0.0 | // ------------------- // | 4.0 5.0 | 6.0 0.0 | // | 7.0 0.0 | 8.0 9.0 | // ------------------- // Number of rows and columns int m = 4; // Number of block rows and block columns int mb = 2; int nb = 2; // BSR block dimension int bsr_dim = 2; // Number of non-zero blocks int nnzb = 3; // BSR row pointers int hbsrRowPtr[3] = {0, 1, 3}; // BSR column indices int hbsrColInd[3] = {0, 0, 1}; // BSR values double hbsrVal[12] = {1.0, 2.0, 0.0, 3.0, 4.0, 7.0, 5.0, 0.0, 6.0, 8.0, 0.0, 9.0}; // Storage scheme of the BSR blocks hipsparseDirection_t dir = HIPSPARSE_DIRECTION_COLUMN; // Transposition of the matrix and rhs matrix hipsparseOperation_t trans = HIPSPARSE_OPERATION_NON_TRANSPOSE; // Solve policy hipsparseSolvePolicy_t solve_policy = HIPSPARSE_SOLVE_POLICY_USE_LEVEL; // Scalar alpha and beta double alpha = 3.7; double hx[4] = {1, 2, 3, 4}; double hy[4]; // Offload data to device int* dbsrRowPtr; int* dbsrColInd; double* dbsrVal; double* dx; double* dy; hipMalloc((void**)&dbsrRowPtr, sizeof(int) * (mb + 1)); hipMalloc((void**)&dbsrColInd, sizeof(int) * nnzb); hipMalloc((void**)&dbsrVal, sizeof(double) * nnzb * bsr_dim * bsr_dim); hipMalloc((void**)&dx, sizeof(double) * nb * bsr_dim); hipMalloc((void**)&dy, sizeof(double) * mb * bsr_dim); hipMemcpy(dbsrRowPtr, hbsrRowPtr, sizeof(int) * (mb + 1), hipMemcpyHostToDevice); hipMemcpy(dbsrColInd, hbsrColInd, sizeof(int) * nnzb, hipMemcpyHostToDevice); hipMemcpy(dbsrVal, hbsrVal, sizeof(double) * nnzb * bsr_dim * bsr_dim, hipMemcpyHostToDevice); hipMemcpy(dx, hx, sizeof(double) * nb * bsr_dim, hipMemcpyHostToDevice); // Matrix descriptor hipsparseMatDescr_t descr; hipsparseCreateMatDescr(&descr); // Matrix fill mode hipsparseSetMatFillMode(descr, HIPSPARSE_FILL_MODE_LOWER); // Matrix diagonal type hipsparseSetMatDiagType(descr, HIPSPARSE_DIAG_TYPE_UNIT); // Matrix info structure bsrsv2Info_t info; hipsparseCreateBsrsv2Info(&info); // Obtain required buffer size int buffer_size; hipsparseDbsrsv2_bufferSize(handle, dir, trans, mb, nnzb, descr, dbsrVal, dbsrRowPtr, dbsrColInd, bsr_dim, info, &buffer_size); // Allocate temporary buffer void* dbuffer; hipMalloc(&dbuffer, buffer_size); // Perform analysis step hipsparseDbsrsv2_analysis(handle, dir, trans, mb, nnzb, descr, dbsrVal, dbsrRowPtr, dbsrColInd, bsr_dim, info, solve_policy, dbuffer); // Call dbsrsm to perform lower triangular solve LX = B hipsparseDbsrsv2_solve(handle, dir, trans, mb, nnzb, &alpha, descr, dbsrVal, dbsrRowPtr, dbsrColInd, bsr_dim, info, dx, dy, solve_policy, dbuffer); // Check for zero pivots int pivot; hipsparseStatus_t status = hipsparseXbsrsv2_zeroPivot(handle, info, &pivot); if(status == HIPSPARSE_STATUS_ZERO_PIVOT) { std::cout << "Found zero pivot in matrix row " << pivot << std::endl; } // Copy results back to the host hipMemcpy(hy, dy, sizeof(double) * mb * bsr_dim, hipMemcpyDeviceToHost); // Clear hipSPARSE hipsparseDestroyBsrsv2Info(info); hipsparseDestroyMatDescr(descr); hipsparseDestroy(handle); // Clear device memory hipFree(dbsrRowPtr); hipFree(dbsrColInd); hipFree(dbsrVal); hipFree(dx); hipFree(dy); hipFree(dbuffer);
Note
The sparse BSR matrix has to be sorted.
Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
Currently, only
trans
== HIPSPARSE_OPERATION_NON_TRANSPOSE andtrans
== HIPSPARSE_OPERATION_TRANSPOSE is supported.
hipsparseXgemvi_bufferSize()#
-
hipsparseStatus_t hipsparseSgemvi_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, int *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseDgemvi_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, int *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseCgemvi_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, int *pBufferSizeInBytes)#
-
hipsparseStatus_t hipsparseZgemvi_bufferSize(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, int nnz, int *pBufferSizeInBytes)#
Dense matrix sparse vector multiplication.
hipsparseXgemvi_bufferSize
returns the size of the temporary storage buffer in bytes required by hipsparseXgemvi(). The temporary storage buffer must be allocated by the user.
hipsparseXgemvi()#
-
hipsparseStatus_t hipsparseSgemvi(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, const float *alpha, const float *A, int lda, int nnz, const float *x, const int *xInd, const float *beta, float *y, hipsparseIndexBase_t idxBase, void *pBuffer)#
-
hipsparseStatus_t hipsparseDgemvi(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, const double *alpha, const double *A, int lda, int nnz, const double *x, const int *xInd, const double *beta, double *y, hipsparseIndexBase_t idxBase, void *pBuffer)#
-
hipsparseStatus_t hipsparseCgemvi(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, const hipComplex *alpha, const hipComplex *A, int lda, int nnz, const hipComplex *x, const int *xInd, const hipComplex *beta, hipComplex *y, hipsparseIndexBase_t idxBase, void *pBuffer)#
-
hipsparseStatus_t hipsparseZgemvi(hipsparseHandle_t handle, hipsparseOperation_t transA, int m, int n, const hipDoubleComplex *alpha, const hipDoubleComplex *A, int lda, int nnz, const hipDoubleComplex *x, const int *xInd, const hipDoubleComplex *beta, hipDoubleComplex *y, hipsparseIndexBase_t idxBase, void *pBuffer)#
Dense matrix sparse vector multiplication.
hipsparseXgemvi
multiplies the scalar \(\alpha\) with a dense \(m \times n\) matrix \(A\) and the sparse vector \(x\) and adds the result to the dense vector \(y\) that is multiplied by the scalar \(\beta\), such that\[ y := \alpha \cdot op(A) \cdot x + \beta \cdot y, \]with\[\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}\]hipsparseXgemvi
requires a user allocated temporary buffer. Its size is returned by hipsparseXgemvi_bufferSize().Note
This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished.
Note
Currently, only
trans
== HIPSPARSE_OPERATION_NON_TRANSPOSE is supported.