This page contains proposed changes for a future release of ROCm. Read the latest Linux release of ROCm documentation for your production environments.

Preconditioner Functions

Contents

Preconditioner Functions#

This module holds all sparse preconditioners.

The sparse preconditioners describe manipulations on a matrix in sparse format to obtain a sparse preconditioner matrix.

hipsparseXbsrilu02_zeroPivot()#

hipsparseStatus_t hipsparseXbsrilu02_zeroPivot(hipsparseHandle_t handle, bsrilu02Info_t info, int *position)#

Incomplete LU factorization with 0 fill-ins and no pivoting using BSR storage format.

hipsparseXbsrilu02_zeroPivot returns HIPSPARSE_STATUS_ZERO_PIVOT, if either a structural or numerical zero has been found during hipsparseXbsrilu02_analysis() or hipsparseXbsrilu02() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored in position, 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

If a zero pivot is found, position \(=j\) means that either the diagonal block \(A_{j,j}\) is missing (structural zero) or the diagonal block \(A_{j,j}\) is not invertible (numerical zero).

Note

hipsparseXbsrilu02_zeroPivot is a blocking function. It might influence performance negatively.

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

  • info[in] structure that holds the information collected during the analysis step.

  • position[inout] pointer to zero pivot \(j\), can be in host or device memory.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, info or position pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_ZERO_PIVOT – zero pivot has been found.

hipsparseXbsrilu02_numericBoost()#

hipsparseStatus_t hipsparseSbsrilu02_numericBoost(hipsparseHandle_t handle, bsrilu02Info_t info, int enable_boost, double *tol, float *boost_val)#
hipsparseStatus_t hipsparseDbsrilu02_numericBoost(hipsparseHandle_t handle, bsrilu02Info_t info, int enable_boost, double *tol, double *boost_val)#
hipsparseStatus_t hipsparseCbsrilu02_numericBoost(hipsparseHandle_t handle, bsrilu02Info_t info, int enable_boost, double *tol, hipComplex *boost_val)#
hipsparseStatus_t hipsparseZbsrilu02_numericBoost(hipsparseHandle_t handle, bsrilu02Info_t info, int enable_boost, double *tol, hipDoubleComplex *boost_val)#

Incomplete LU factorization with 0 fill-ins and no pivoting using BSR storage format.

hipsparseXbsrilu02_numericBoost enables the user to replace a numerical value in an incomplete LU factorization. tol is used to determine whether a numerical value is replaced by boost_val, such that \(A_{j,j} = \text{boost_val}\) if \(\text{tol} \ge \left|A_{j,j}\right|\).

Note

The boost value is enabled by setting enable_boost to 1 or disabled by setting enable_boost to 0.

Note

tol and boost_val can be in host or device memory.

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

  • info[in] structure that holds the information collected during the analysis step.

  • enable_boost[in] enable/disable numeric boost.

  • tol[in] tolerance to determine whether a numerical value is replaced or not.

  • boost_val[in] boost value to replace a numerical value.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, info, tol or boost_val pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXbsrilu02_bufferSize()#

hipsparseStatus_t hipsparseSbsrilu02_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, float *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrilu02Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDbsrilu02_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, double *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrilu02Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCbsrilu02_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrilu02Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZbsrilu02_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipDoubleComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrilu02Info_t info, int *pBufferSizeInBytes)#

Incomplete LU factorization with 0 fill-ins and no pivoting using BSR storage format.

hipsparseXbsrilu02_bufferSize returns the size of the temporary storage buffer in bytes that is required by hipsparseXbsrilu02_analysis() and hipsparseXbsrilu02(). The temporary storage buffer must be allocated by the user.

Parameters:
Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, mb, nnzb, blockDim, descrA, bsrSortedValA, bsrSortedRowPtrA, bsrSortedColIndA, info or pBufferSizeInBytes pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXbsrilu02_analysis()#

hipsparseStatus_t hipsparseSbsrilu02_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, float *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseDbsrilu02_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, double *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseCbsrilu02_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseZbsrilu02_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipDoubleComplex *bsrSortedValA, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#

Incomplete LU factorization with 0 fill-ins and no pivoting using BSR storage format.

hipsparseXbsrilu02_analysis performs the analysis step for hipsparseXbsrilu02().

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.

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

  • dirA[in] direction that specified whether to count nonzero elements by HIPSPARSE_DIRECTION_ROW or by HIPSPARSE_DIRECTION_COLUMN.

  • mb[in] number of block rows in the sparse BSR matrix.

  • nnzb[in] number of non-zero block entries of the sparse BSR matrix.

  • descrA[in] descriptor of the sparse BSR matrix.

  • bsrSortedValA[in] array of length nnzb*blockDim*blockDim containing the values of the sparse BSR matrix.

  • bsrSortedRowPtrA[in] array of mb+1 elements that point to the start of every block row of the sparse BSR matrix.

  • bsrSortedColIndA[in] array of nnzb elements containing the block column indices of the sparse BSR matrix.

  • blockDim[in] the block dimension of the BSR matrix. Between 1 and m where m=mb*blockDim.

  • info[out] structure that holds the information collected during the analysis step.

  • policy[in] HIPSPARSE_SOLVE_POLICY_NO_LEVEL or HIPSPARSE_SOLVE_POLICY_USE_LEVEL.

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, mb, nnzb, blockDim, descrA, bsrSortedValA, bsrSortedRowPtrA, bsrSortedColIndA, info or pBuffer pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXbsrilu02()#

hipsparseStatus_t hipsparseSbsrilu02(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, float *bsrSortedValA_valM, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseDbsrilu02(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, double *bsrSortedValA_valM, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseCbsrilu02(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipComplex *bsrSortedValA_valM, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseZbsrilu02(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipDoubleComplex *bsrSortedValA_valM, const int *bsrSortedRowPtrA, const int *bsrSortedColIndA, int blockDim, bsrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#

Incomplete LU factorization with 0 fill-ins and no pivoting using BSR storage format.

hipsparseXbsrilu02 computes the incomplete LU factorization with 0 fill-ins and no pivoting of a sparse \(mb \times mb\) BSR matrix \(A\), such that

\[ A \approx LU \]

hipsparseXbsrilu02 requires a user allocated temporary buffer. Its size is returned by hipsparseXbsrilu02_bufferSize(). Furthermore, analysis meta data is required. It can be obtained by hipsparseXbsrilu02_analysis(). hipsparseXbsrilu02 reports the first zero pivot (either numerical or structural zero). The zero pivot status can be obtained by calling hipsparseXbsrilu02_zeroPivot().

Note

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

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

  • dirA[in] direction that specified whether to count nonzero elements by HIPSPARSE_DIRECTION_ROW or by HIPSPARSE_DIRECTION_COLUMN.

  • mb[in] number of block rows in the sparse BSR matrix.

  • nnzb[in] number of non-zero block entries of the sparse BSR matrix.

  • descrA[in] descriptor of the sparse BSR matrix.

  • bsrSortedValA_valM[inout] array of length nnzb*blockDim*blockDim containing the values of the sparse BSR matrix.

  • bsrSortedRowPtrA[in] array of mb+1 elements that point to the start of every block row of the sparse BSR matrix.

  • bsrSortedColIndA[in] array of nnzb elements containing the block column indices of the sparse BSR matrix.

  • blockDim[in] the block dimension of the BSR matrix. Between 1 and m where m=mb*blockDim.

  • info[in] structure that holds the information collected during the analysis step.

  • policy[in] HIPSPARSE_SOLVE_POLICY_NO_LEVEL or HIPSPARSE_SOLVE_POLICY_USE_LEVEL.

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, mb, nnzb, blockDim, descrA, bsrSortedValA_valM, bsrSortedRowPtrA or bsrSortedColIndA pointer is invalid.

  • HIPSPARSE_STATUS_ARCH_MISMATCH – the device is not supported.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXcsrilu02_zeroPivot()#

hipsparseStatus_t hipsparseXcsrilu02_zeroPivot(hipsparseHandle_t handle, csrilu02Info_t info, int *position)#

Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format.

hipsparseXcsrilu02_zeroPivot returns HIPSPARSE_STATUS_ZERO_PIVOT, if either a structural or numerical zero has been found during hipsparseXcsrilu02() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored in position, 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

hipsparseXcsrilu02_zeroPivot is a blocking function. It might influence performance negatively.

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

  • info[in] structure that holds the information collected during the analysis step.

  • position[inout] pointer to zero pivot \(j\), can be in host or device memory.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, info or position pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_ZERO_PIVOT – zero pivot has been found.

hipsparseXcsrilu02_numericBoost()#

hipsparseStatus_t hipsparseScsrilu02_numericBoost(hipsparseHandle_t handle, csrilu02Info_t info, int enable_boost, double *tol, float *boost_val)#
hipsparseStatus_t hipsparseDcsrilu02_numericBoost(hipsparseHandle_t handle, csrilu02Info_t info, int enable_boost, double *tol, double *boost_val)#
hipsparseStatus_t hipsparseCcsrilu02_numericBoost(hipsparseHandle_t handle, csrilu02Info_t info, int enable_boost, double *tol, hipComplex *boost_val)#
hipsparseStatus_t hipsparseZcsrilu02_numericBoost(hipsparseHandle_t handle, csrilu02Info_t info, int enable_boost, double *tol, hipDoubleComplex *boost_val)#

Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format.

hipsparseXcsrilu02_numericBoost enables the user to replace a numerical value in an incomplete LU factorization. tol is used to determine whether a numerical value is replaced by boost_val, such that \(A_{j,j} = \text{boost_val}\) if \(\text{tol} \ge \left|A_{j,j}\right|\).

Note

The boost value is enabled by setting enable_boost to 1 or disabled by setting enable_boost to 0.

Note

tol and boost_val can be in host or device memory.

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

  • info[in] structure that holds the information collected during the analysis step.

  • enable_boost[in] enable/disable numeric boost.

  • tol[in] tolerance to determine whether a numerical value is replaced or not.

  • boost_val[in] boost value to replace a numerical value.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, info, tol or boost_val pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXcsrilu02_bufferSize()#

hipsparseStatus_t hipsparseScsrilu02_bufferSize(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDcsrilu02_bufferSize(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCcsrilu02_bufferSize(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZcsrilu02_bufferSize(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, int *pBufferSizeInBytes)#

Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format.

hipsparseXcsrilu02_bufferSize returns the size of the temporary storage buffer in bytes that is required by hipsparseXcsrilu02_analysis() and hipsparseXcsrilu02(). The temporary storage buffer must be allocated by the user.

Parameters:
Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, nnz, descrA, csrSortedValA, csrSortedRowPtrA, csrSortedColIndA, info or pBufferSizeInBytes pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXcsrilu02_bufferSizeExt()#

hipsparseStatus_t hipsparseScsrilu02_bufferSizeExt(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDcsrilu02_bufferSizeExt(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCcsrilu02_bufferSizeExt(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZcsrilu02_bufferSizeExt(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, size_t *pBufferSizeInBytes)#

Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format.

hipsparseXcsrilu02_bufferSizeExt returns the size of the temporary storage buffer in bytes that is required by hipsparseXcsrilu02_analysis() and hipsparseXcsrilu02(). The temporary storage buffer must be allocated by the user.

Parameters:
Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, nnz, descrA, csrSortedValA, csrSortedRowPtrA, csrSortedColIndA, info or pBufferSizeInBytes pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXcsrilu02_analysis()#

hipsparseStatus_t hipsparseScsrilu02_analysis(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, const float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseDcsrilu02_analysis(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, const double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseCcsrilu02_analysis(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, const hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseZcsrilu02_analysis(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, const hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#

Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format.

hipsparseXcsrilu02_analysis performs the analysis step for hipsparseXcsrilu02().

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.

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

  • m[in] number of rows of the sparse CSR matrix.

  • nnz[in] number of non-zero entries of the sparse CSR matrix.

  • descrA[in] descriptor of the sparse CSR matrix.

  • csrSortedValA[in] array of nnz elements of the sparse CSR matrix.

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

  • csrSortedColIndA[in] array of nnz elements containing the column indices of the sparse CSR matrix.

  • info[out] structure that holds the information collected during the analysis step.

  • policy[in] HIPSPARSE_SOLVE_POLICY_NO_LEVEL or HIPSPARSE_SOLVE_POLICY_USE_LEVEL.

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, nnz, descrA, csrSortedValA, csrSortedRowPtrA, csrSortedColIndA, info or pBuffer pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXcsrilu02()#

hipsparseStatus_t hipsparseScsrilu02(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, float *csrSortedValA_valM, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseDcsrilu02(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, double *csrSortedValA_valM, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseCcsrilu02(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, hipComplex *csrSortedValA_valM, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseZcsrilu02(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, hipDoubleComplex *csrSortedValA_valM, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csrilu02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#

Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format.

hipsparseXcsrilu02 computes the incomplete LU factorization with 0 fill-ins and no pivoting of a sparse \(m \times m\) CSR matrix \(A\), such that

\[ A \approx LU \]

hipsparseXcsrilu02 requires a user allocated temporary buffer. Its size is returned by hipsparseXcsrilu02_bufferSize() or hipsparseXcsrilu02_bufferSizeExt(). Furthermore, analysis meta data is required. It can be obtained by hipsparseXcsrilu02_analysis(). hipsparseXcsrilu02 reports the first zero pivot (either numerical or structural zero). The zero pivot status can be obtained by calling hipsparseXcsrilu02_zeroPivot().

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.

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

  • m[in] number of rows of the sparse CSR matrix.

  • nnz[in] number of non-zero entries of the sparse CSR matrix.

  • descrA[in] descriptor of the sparse CSR matrix.

  • csrSortedValA_valM[inout] array of nnz elements of the sparse CSR matrix.

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

  • csrSortedColIndA[in] array of nnz elements containing the column indices of the sparse CSR matrix.

  • info[in] structure that holds the information collected during the analysis step.

  • policy[in] HIPSPARSE_SOLVE_POLICY_NO_LEVEL or HIPSPARSE_SOLVE_POLICY_USE_LEVEL.

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, nnz, descrA, csrSortedValA_valM, csrSortedRowPtrA or csrSortedColIndA pointer is invalid.

  • HIPSPARSE_STATUS_ARCH_MISMATCH – the device is not supported.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXbsric02_zeroPivot()#

hipsparseStatus_t hipsparseXbsric02_zeroPivot(hipsparseHandle_t handle, bsric02Info_t info, int *position)#

Incomplete Cholesky factorization with 0 fill-ins and no pivoting using BSR storage format.

hipsparseXbsric02_zeroPivot returns HIPSPARSE_STATUS_ZERO_PIVOT, if either a structural or numerical zero has been found during hipsparseXbsric02_analysis() or hipsparseXbsric02() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored in position, 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

If a zero pivot is found, position=j means that either the diagonal block A(j,j) is missing (structural zero) or the diagonal block A(j,j) is not positive definite (numerical zero).

Note

hipsparseXbsric02_zeroPivot is a blocking function. It might influence performance negatively.

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

  • info[in] structure that holds the information collected during the analysis step.

  • position[inout] pointer to zero pivot \(j\), can be in host or device memory.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, info or position pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_ZERO_PIVOT – zero pivot has been found.

hipsparseXbsric02_bufferSize()#

hipsparseStatus_t hipsparseSbsric02_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, float *bsrValA, const int *bsrRowPtrA, const int *bsrColIndA, int blockDim, bsric02Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDbsric02_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, double *bsrValA, const int *bsrRowPtrA, const int *bsrColIndA, int blockDim, bsric02Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCbsric02_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipComplex *bsrValA, const int *bsrRowPtrA, const int *bsrColIndA, int blockDim, bsric02Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZbsric02_bufferSize(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipDoubleComplex *bsrValA, const int *bsrRowPtrA, const int *bsrColIndA, int blockDim, bsric02Info_t info, int *pBufferSizeInBytes)#

Incomplete Cholesky factorization with 0 fill-ins and no pivoting using BSR storage format.

hipsparseXbsric02_bufferSize returns the size of the temporary storage buffer in bytes that is required by hipsparseXbsric02_analysis() and hipsparseXbsric02(). The temporary storage buffer must be allocated by the user.

Parameters:
Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, mb, nnzb, blockDim, descrA, bsrValA, bsrRowPtrA, bsrColIndA, info or pBufferSizeInBytes pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXbsric02_analysis()#

hipsparseStatus_t hipsparseSbsric02_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, const float *bsrValA, const int *bsrRowPtrA, const int *bsrColIndA, int blockDim, bsric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseDbsric02_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, const double *bsrValA, const int *bsrRowPtrA, const int *bsrColIndA, int blockDim, bsric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseCbsric02_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, const hipComplex *bsrValA, const int *bsrRowPtrA, const int *bsrColIndA, int blockDim, bsric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseZbsric02_analysis(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, const hipDoubleComplex *bsrValA, const int *bsrRowPtrA, const int *bsrColIndA, int blockDim, bsric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#

Incomplete Cholesky factorization with 0 fill-ins and no pivoting using BSR storage format.

hipsparseXbsric02_analysis performs the analysis step for hipsparseXbsric02().

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.

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

  • dirA[in] direction that specified whether to count nonzero elements by HIPSPARSE_DIRECTION_ROW or by HIPSPARSE_DIRECTION_COLUMN.

  • mb[in] number of block rows in the sparse BSR matrix.

  • nnzb[in] number of non-zero block entries of the sparse BSR matrix.

  • descrA[in] descriptor of the sparse BSR matrix.

  • bsrValA[in] array of length nnzb*blockDim*blockDim containing the values of the sparse BSR matrix.

  • bsrRowPtrA[in] array of mb+1 elements that point to the start of every block row of the sparse BSR matrix.

  • bsrColIndA[in] array of nnzb elements containing the block column indices of the sparse BSR matrix.

  • blockDim[in] the block dimension of the BSR matrix. Between 1 and m where m=mb*blockDim.

  • info[out] structure that holds the information collected during the analysis step.

  • policy[in] HIPSPARSE_SOLVE_POLICY_NO_LEVEL or HIPSPARSE_SOLVE_POLICY_USE_LEVEL.

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, mb, nnzb, blockDim, descrA, bsrValA, bsrRowPtrA, bsrColIndA, info or pBuffer pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXbsric02()#

hipsparseStatus_t hipsparseSbsric02(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, float *bsrValA, const int *bsrRowPtrA, const int *bsrColIndA, int blockDim, bsric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseDbsric02(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, double *bsrValA, const int *bsrRowPtrA, const int *bsrColIndA, int blockDim, bsric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseCbsric02(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipComplex *bsrValA, const int *bsrRowPtrA, const int *bsrColIndA, int blockDim, bsric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseZbsric02(hipsparseHandle_t handle, hipsparseDirection_t dirA, int mb, int nnzb, const hipsparseMatDescr_t descrA, hipDoubleComplex *bsrValA, const int *bsrRowPtrA, const int *bsrColIndA, int blockDim, bsric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#

Incomplete Cholesky factorization with 0 fill-ins and no pivoting using BSR storage format.

hipsparseXbsric02 computes the incomplete Cholesky factorization with 0 fill-ins and no pivoting of a sparse \(mb \times mb\) BSR matrix \(A\), such that

\[ A \approx LL^T \]

hipsparseXbsric02 requires a user allocated temporary buffer. Its size is returned by hipsparseXbsric02_bufferSize(). Furthermore, analysis meta data is required. It can be obtained by hipsparseXbsric02_analysis(). hipsparseXbsric02 reports the first zero pivot (either numerical or structural zero). The zero pivot status can be obtained by calling hipsparseXbsric02_zeroPivot().

Note

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

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

  • dirA[in] direction that specified whether to count nonzero elements by HIPSPARSE_DIRECTION_ROW or by HIPSPARSE_DIRECTION_COLUMN.

  • mb[in] number of block rows in the sparse BSR matrix.

  • nnzb[in] number of non-zero block entries of the sparse BSR matrix.

  • descrA[in] descriptor of the sparse BSR matrix.

  • bsrValA[inout] array of length nnzb*blockDim*blockDim containing the values of the sparse BSR matrix.

  • bsrRowPtrA[in] array of mb+1 elements that point to the start of every block row of the sparse BSR matrix.

  • bsrColIndA[in] array of nnzb elements containing the block column indices of the sparse BSR matrix.

  • blockDim[in] the block dimension of the BSR matrix. Between 1 and m where m=mb*blockDim.

  • info[in] structure that holds the information collected during the analysis step.

  • policy[in] HIPSPARSE_SOLVE_POLICY_NO_LEVEL or HIPSPARSE_SOLVE_POLICY_USE_LEVEL.

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, mb, nnzb, blockDim, descrA, bsrValA, bsrRowPtrA, or bsrColIndA pointer is invalid.

  • HIPSPARSE_STATUS_ARCH_MISMATCH – the device is not supported.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXcsric02_zeroPivot()#

hipsparseStatus_t hipsparseXcsric02_zeroPivot(hipsparseHandle_t handle, csric02Info_t info, int *position)#

Incomplete Cholesky factorization with 0 fill-ins and no pivoting using CSR storage format.

hipsparseXcsric02_zeroPivot returns HIPSPARSE_STATUS_ZERO_PIVOT, if either a structural or numerical zero has been found during hipsparseXcsric02_analysis() or hipsparseXcsric02() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored in position, 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

hipsparseXcsric02_zeroPivot is a blocking function. It might influence performance negatively.

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

  • info[in] structure that holds the information collected during the analysis step.

  • position[inout] pointer to zero pivot \(j\), can be in host or device memory.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, info or position pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_ZERO_PIVOT – zero pivot has been found.

hipsparseXcsric02_bufferSize()#

hipsparseStatus_t hipsparseScsric02_bufferSize(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDcsric02_bufferSize(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCcsric02_bufferSize(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, int *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZcsric02_bufferSize(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, int *pBufferSizeInBytes)#

Incomplete Cholesky factorization with 0 fill-ins and no pivoting using CSR storage format.

hipsparseXcsric02_bufferSize returns the size of the temporary storage buffer in bytes that is required by hipsparseXcsric02_analysis() and hipsparseXcsric02().

Parameters:
Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, nnz, descrA, csrSortedValA, csrSortedRowPtrA, csrSortedColIndA, info or pBufferSizeInBytes pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXcsric02_bufferSizeExt()#

hipsparseStatus_t hipsparseScsric02_bufferSizeExt(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDcsric02_bufferSizeExt(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCcsric02_bufferSizeExt(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZcsric02_bufferSizeExt(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, size_t *pBufferSizeInBytes)#

Incomplete Cholesky factorization with 0 fill-ins and no pivoting using CSR storage format.

hipsparseXcsric02_bufferSizeExt returns the size of the temporary storage buffer in bytes that is required by hipsparseXcsric02_analysis() and hipsparseXcsric02().

Parameters:
Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, nnz, descrA, csrSortedValA, csrSortedRowPtrA, csrSortedColIndA, info or pBufferSizeInBytes pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXcsric02_analysis()#

hipsparseStatus_t hipsparseScsric02_analysis(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, const float *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseDcsric02_analysis(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, const double *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseCcsric02_analysis(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, const hipComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseZcsric02_analysis(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, const hipDoubleComplex *csrSortedValA, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#

Incomplete Cholesky factorization with 0 fill-ins and no pivoting using CSR storage format.

hipsparseXcsric02_analysis performs the analysis step for hipsparseXcsric02().

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.

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

  • m[in] number of rows of the sparse CSR matrix.

  • nnz[in] number of non-zero entries of the sparse CSR matrix.

  • descrA[in] descriptor of the sparse CSR matrix.

  • csrSortedValA[in] array of nnz elements of the sparse CSR matrix.

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

  • csrSortedColIndA[in] array of nnz elements containing the column indices of the sparse CSR matrix.

  • info[out] structure that holds the information collected during the analysis step.

  • policy[in] HIPSPARSE_SOLVE_POLICY_NO_LEVEL or HIPSPARSE_SOLVE_POLICY_USE_LEVEL.

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, nnz, descrA, csrSortedValA, csrSortedRowPtrA, csrSortedColIndA, info or pBuffer pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXcsric02()#

hipsparseStatus_t hipsparseScsric02(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, float *csrSortedValA_valM, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseDcsric02(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, double *csrSortedValA_valM, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseCcsric02(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, hipComplex *csrSortedValA_valM, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#
hipsparseStatus_t hipsparseZcsric02(hipsparseHandle_t handle, int m, int nnz, const hipsparseMatDescr_t descrA, hipDoubleComplex *csrSortedValA_valM, const int *csrSortedRowPtrA, const int *csrSortedColIndA, csric02Info_t info, hipsparseSolvePolicy_t policy, void *pBuffer)#

Incomplete Cholesky factorization with 0 fill-ins and no pivoting using CSR storage format.

hipsparseXcsric02 computes the incomplete Cholesky factorization with 0 fill-ins and no pivoting of a sparse \(m \times m\) CSR matrix \(A\), such that

\[ A \approx LL^T \]

hipsparseXcsric02 requires a user allocated temporary buffer. Its size is returned by hipsparseXcsric02_bufferSize() or hipsparseXcsric02_bufferSizeExt(). Furthermore, analysis meta data is required. It can be obtained by hipsparseXcsric02_analysis(). hipsparseXcsric02 reports the first zero pivot (either numerical or structural zero). The zero pivot status can be obtained by calling hipsparseXcsric02_zeroPivot().

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.

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

  • m[in] number of rows of the sparse CSR matrix.

  • nnz[in] number of non-zero entries of the sparse CSR matrix.

  • descrA[in] descriptor of the sparse CSR matrix.

  • csrSortedValA_valM[inout] array of nnz elements of the sparse CSR matrix.

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

  • csrSortedColIndA[in] array of nnz elements containing the column indices of the sparse CSR matrix.

  • info[in] structure that holds the information collected during the analysis step.

  • policy[in] HIPSPARSE_SOLVE_POLICY_NO_LEVEL or HIPSPARSE_SOLVE_POLICY_USE_LEVEL.

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, nnz, descrA, csrSortedValA_valM, csrSortedRowPtrA or csrSortedColIndA pointer is invalid.

  • HIPSPARSE_STATUS_ARCH_MISMATCH – the device is not supported.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

  • HIPSPARSE_STATUS_NOT_SUPPORTEDhipsparseMatrixType_t != HIPSPARSE_MATRIX_TYPE_GENERAL.

hipsparseXgtsv2_bufferSizeExt()#

hipsparseStatus_t hipsparseSgtsv2_bufferSizeExt(hipsparseHandle_t handle, int m, int n, const float *dl, const float *d, const float *du, const float *B, int ldb, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDgtsv2_bufferSizeExt(hipsparseHandle_t handle, int m, int n, const double *dl, const double *d, const double *du, const double *B, int db, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCgtsv2_bufferSizeExt(hipsparseHandle_t handle, int m, int n, const hipComplex *dl, const hipComplex *d, const hipComplex *du, const hipComplex *B, int ldb, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZgtsv2_bufferSizeExt(hipsparseHandle_t handle, int m, int n, const hipDoubleComplex *dl, const hipDoubleComplex *d, const hipDoubleComplex *du, const hipDoubleComplex *B, int ldb, size_t *pBufferSizeInBytes)#

Tridiagonal solver with pivoting.

hipsparseXgtsv2_bufferSize returns the size of the temporary storage buffer in bytes that is required by hipsparseXgtsv2(). The temporary storage buffer must be allocated by the user.

hipsparseXgtsv2()#

hipsparseStatus_t hipsparseSgtsv2(hipsparseHandle_t handle, int m, int n, const float *dl, const float *d, const float *du, float *B, int ldb, void *pBuffer)#
hipsparseStatus_t hipsparseDgtsv2(hipsparseHandle_t handle, int m, int n, const double *dl, const double *d, const double *du, double *B, int ldb, void *pBuffer)#
hipsparseStatus_t hipsparseCgtsv2(hipsparseHandle_t handle, int m, int n, const hipComplex *dl, const hipComplex *d, const hipComplex *du, hipComplex *B, int ldb, void *pBuffer)#
hipsparseStatus_t hipsparseZgtsv2(hipsparseHandle_t handle, int m, int n, const hipDoubleComplex *dl, const hipDoubleComplex *d, const hipDoubleComplex *du, hipDoubleComplex *B, int ldb, void *pBuffer)#

Tridiagonal solver with pivoting.

hipsparseXgtsv2 solves a tridiagonal system for multiple right hand sides using pivoting.

Note

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

hipsparseXgtsv2_nopivot_bufferSizeExt()#

hipsparseStatus_t hipsparseSgtsv2_nopivot_bufferSizeExt(hipsparseHandle_t handle, int m, int n, const float *dl, const float *d, const float *du, const float *B, int ldb, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDgtsv2_nopivot_bufferSizeExt(hipsparseHandle_t handle, int m, int n, const double *dl, const double *d, const double *du, const double *B, int ldb, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCgtsv2_nopivot_bufferSizeExt(hipsparseHandle_t handle, int m, int n, const hipComplex *dl, const hipComplex *d, const hipComplex *du, const hipComplex *B, int ldb, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZgtsv2_nopivot_bufferSizeExt(hipsparseHandle_t handle, int m, int n, const hipDoubleComplex *dl, const hipDoubleComplex *d, const hipDoubleComplex *du, const hipDoubleComplex *B, int ldb, size_t *pBufferSizeInBytes)#

Tridiagonal solver (no pivoting)

hipsparseXgtsv2_nopivot_bufferSizeExt returns the size of the temporary storage buffer in bytes that is required by hipsparseXgtsv2_nopivot(). The temporary storage buffer must be allocated by the user.

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

  • m[in] size of the tri-diagonal linear system (must be >= 2).

  • n[in] number of columns in the dense matrix B.

  • dl[in] lower diagonal of tri-diagonal system. First entry must be zero.

  • d[in] main diagonal of tri-diagonal system.

  • du[in] upper diagonal of tri-diagonal system. Last entry must be zero.

  • B[in] Dense matrix of size ( ldb, n ).

  • ldb[in] Leading dimension of B. Must satisfy ldb >= max(1, m).

  • pBufferSizeInBytes[out] number of bytes of the temporary storage buffer required by hipsparseSgtsv2_nopivot(), hipsparseDgtsv2_nopivot(), hipsparseCgtsv2_nopivot() and hipsparseZgtsv2_nopivot().

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, n, ldb, dl, d, du, B or pBufferSizeInBytes pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXgtsv2_nopivot()#

hipsparseStatus_t hipsparseSgtsv2_nopivot(hipsparseHandle_t handle, int m, int n, const float *dl, const float *d, const float *du, float *B, int ldb, void *pBuffer)#
hipsparseStatus_t hipsparseDgtsv2_nopivot(hipsparseHandle_t handle, int m, int n, const double *dl, const double *d, const double *du, double *B, int ldb, void *pBuffer)#
hipsparseStatus_t hipsparseCgtsv2_nopivot(hipsparseHandle_t handle, int m, int n, const hipComplex *dl, const hipComplex *d, const hipComplex *du, hipComplex *B, int ldb, void *pBuffer)#
hipsparseStatus_t hipsparseZgtsv2_nopivot(hipsparseHandle_t handle, int m, int n, const hipDoubleComplex *dl, const hipDoubleComplex *d, const hipDoubleComplex *du, hipDoubleComplex *B, int ldb, void *pBuffer)#

Tridiagonal solver (no pivoting)

hipsparseXgtsv2_nopivot solves a tridiagonal linear system for multiple right-hand sides

Note

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

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

  • m[in] size of the tri-diagonal linear system (must be >= 2).

  • n[in] number of columns in the dense matrix B.

  • dl[in] lower diagonal of tri-diagonal system. First entry must be zero.

  • d[in] main diagonal of tri-diagonal system.

  • du[in] upper diagonal of tri-diagonal system. Last entry must be zero.

  • B[inout] Dense matrix of size ( ldb, n ).

  • ldb[in] Leading dimension of B. Must satisfy ldb >= max(1, m).

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, n, ldb, dl, d, du, B or pBuffer pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXgtsv2StridedBatch_bufferSizeExt()#

hipsparseStatus_t hipsparseSgtsv2StridedBatch_bufferSizeExt(hipsparseHandle_t handle, int m, const float *dl, const float *d, const float *du, const float *x, int batchCount, int batchStride, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDgtsv2StridedBatch_bufferSizeExt(hipsparseHandle_t handle, int m, const double *dl, const double *d, const double *du, const double *x, int batchCount, int batchStride, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCgtsv2StridedBatch_bufferSizeExt(hipsparseHandle_t handle, int m, const hipComplex *dl, const hipComplex *d, const hipComplex *du, const hipComplex *x, int batchCount, int batchStride, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZgtsv2StridedBatch_bufferSizeExt(hipsparseHandle_t handle, int m, const hipDoubleComplex *dl, const hipDoubleComplex *d, const hipDoubleComplex *du, const hipDoubleComplex *x, int batchCount, int batchStride, size_t *pBufferSizeInBytes)#

Strided Batch tridiagonal solver (no pivoting)

hipsparseXgtsv2StridedBatch_bufferSizeExt returns the size of the temporary storage buffer in bytes that is required by hipsparseXgtsv2StridedBatch(). The temporary storage buffer must be allocated by the user.

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

  • m[in] size of the tri-diagonal linear system.

  • dl[in] lower diagonal of tri-diagonal system where the ith system lower diagonal starts at dl+batchStride*i.

  • d[in] main diagonal of tri-diagonal system where the ith system diagonal starts at d+batchStride*i.

  • du[in] upper diagonal of tri-diagonal system where the ith system upper diagonal starts at du+batchStride*i.

  • x[inout] Dense array of righthand-sides where the ith righthand-side starts at x+batchStride*i.

  • batchCount[in] The number of systems to solve.

  • batchStride[in] The number of elements that separate each system. Must satisfy batchStride >= m.

  • pBufferSizeInBytes[out] number of bytes of the temporary storage buffer required by hipsparseSgtsv2StridedBatch(), hipsparseDgtsv2StridedBatch(), hipsparseCgtsv2StridedBatch() and hipsparseZgtsv2StridedBatch().

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, batchCount, batchStride, dl, d, du, x or pBufferSizeInBytes pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXgtsv2StridedBatch()#

hipsparseStatus_t hipsparseSgtsv2StridedBatch(hipsparseHandle_t handle, int m, const float *dl, const float *d, const float *du, float *x, int batchCount, int batchStride, void *pBuffer)#
hipsparseStatus_t hipsparseDgtsv2StridedBatch(hipsparseHandle_t handle, int m, const double *dl, const double *d, const double *du, double *x, int batchCount, int batchStride, void *pBuffer)#
hipsparseStatus_t hipsparseCgtsv2StridedBatch(hipsparseHandle_t handle, int m, const hipComplex *dl, const hipComplex *d, const hipComplex *du, hipComplex *x, int batchCount, int batchStride, void *pBuffer)#
hipsparseStatus_t hipsparseZgtsv2StridedBatch(hipsparseHandle_t handle, int m, const hipDoubleComplex *dl, const hipDoubleComplex *d, const hipDoubleComplex *du, hipDoubleComplex *x, int batchCount, int batchStride, void *pBuffer)#

Strided Batch tridiagonal solver (no pivoting)

hipsparseXgtsv2StridedBatch solves a batched tridiagonal linear system

Note

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

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

  • m[in] size of the tri-diagonal linear system (must be >= 2).

  • dl[in] lower diagonal of tri-diagonal system. First entry must be zero.

  • d[in] main diagonal of tri-diagonal system.

  • du[in] upper diagonal of tri-diagonal system. Last entry must be zero.

  • x[inout] Dense array of righthand-sides where the ith righthand-side starts at x+batchStride*i.

  • batchCount[in] The number of systems to solve.

  • batchStride[in] The number of elements that separate each system. Must satisfy batchStride >= m.

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, batchCount, batchStride, dl, d, du, x or pBuffer pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXgtsvInterleavedBatch_bufferSizeExt()#

hipsparseStatus_t hipsparseSgtsvInterleavedBatch_bufferSizeExt(hipsparseHandle_t handle, int algo, int m, const float *dl, const float *d, const float *du, const float *x, int batchCount, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDgtsvInterleavedBatch_bufferSizeExt(hipsparseHandle_t handle, int algo, int m, const double *dl, const double *d, const double *du, const double *x, int batchCount, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCgtsvInterleavedBatch_bufferSizeExt(hipsparseHandle_t handle, int algo, int m, const hipComplex *dl, const hipComplex *d, const hipComplex *du, const hipComplex *x, int batchCount, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZgtsvInterleavedBatch_bufferSizeExt(hipsparseHandle_t handle, int algo, int m, const hipDoubleComplex *dl, const hipDoubleComplex *d, const hipDoubleComplex *du, const hipDoubleComplex *x, int batchCount, size_t *pBufferSizeInBytes)#

Interleaved Batch tridiagonal solver.

hipsparseXgtsvInterleavedBatch_bufferSizeExt returns the size of the temporary storage buffer in bytes that is required by hipsparseXgtsvInterleavedBatch(). The temporary storage buffer must be allocated by the user.

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

  • algo[in] Algorithm to use when solving tridiagonal systems. Options are thomas ( algo=0 ), LU ( algo=1 ), or QR ( algo=2 ). Thomas algorithm is the fastest but is not stable while LU and QR are slower but are stable.

  • m[in] size of the tri-diagonal linear system.

  • dl[in] lower diagonal of tri-diagonal system. The first element of the lower diagonal must be zero.

  • d[in] main diagonal of tri-diagonal system.

  • du[in] upper diagonal of tri-diagonal system. The last element of the upper diagonal must be zero.

  • x[inout] Dense array of righthand-sides with dimension batchCount by m.

  • batchCount[in] The number of systems to solve.

  • pBufferSizeInBytes[out] number of bytes of the temporary storage buffer required by hipsparseSgtsvInterleavedBatch(), hipsparseDgtsvInterleavedBatch(), hipsparseCgtsvInterleavedBatch() and hipsparseZgtsvInterleavedBatch().

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, batchCount, dl, d, du, x or pBufferSizeInBytes pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXgtsvInterleavedBatch()#

hipsparseStatus_t hipsparseSgtsvInterleavedBatch(hipsparseHandle_t handle, int algo, int m, float *dl, float *d, float *du, float *x, int batchCount, void *pBuffer)#
hipsparseStatus_t hipsparseDgtsvInterleavedBatch(hipsparseHandle_t handle, int algo, int m, double *dl, double *d, double *du, double *x, int batchCount, void *pBuffer)#
hipsparseStatus_t hipsparseCgtsvInterleavedBatch(hipsparseHandle_t handle, int algo, int m, hipComplex *dl, hipComplex *d, hipComplex *du, hipComplex *x, int batchCount, void *pBuffer)#
hipsparseStatus_t hipsparseZgtsvInterleavedBatch(hipsparseHandle_t handle, int algo, int m, hipDoubleComplex *dl, hipDoubleComplex *d, hipDoubleComplex *du, hipDoubleComplex *x, int batchCount, void *pBuffer)#

Interleaved Batch tridiagonal solver.

hipsparseXgtsvInterleavedBatch solves a batched tridiagonal linear system

Note

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

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

  • algo[in] Algorithm to use when solving tridiagonal systems. Options are thomas ( algo=0 ), LU ( algo=1 ), or QR ( algo=2 ). Thomas algorithm is the fastest but is not stable while LU and QR are slower but are stable.

  • m[in] size of the tri-diagonal linear system.

  • dl[inout] lower diagonal of tri-diagonal system. The first element of the lower diagonal must be zero.

  • d[inout] main diagonal of tri-diagonal system.

  • du[inout] upper diagonal of tri-diagonal system. The last element of the upper diagonal must be zero.

  • x[inout] Dense array of righthand-sides with dimension batchCount by m.

  • batchCount[in] The number of systems to solve.

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, batchCount, dl, d, du, x or pBuffer pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXgpsvInterleavedBatch_bufferSizeExt()#

hipsparseStatus_t hipsparseSgpsvInterleavedBatch_bufferSizeExt(hipsparseHandle_t handle, int algo, int m, const float *ds, const float *dl, const float *d, const float *du, const float *dw, const float *x, int batchCount, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseDgpsvInterleavedBatch_bufferSizeExt(hipsparseHandle_t handle, int algo, int m, const double *ds, const double *dl, const double *d, const double *du, const double *dw, const double *x, int batchCount, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseCgpsvInterleavedBatch_bufferSizeExt(hipsparseHandle_t handle, int algo, int m, const hipComplex *ds, const hipComplex *dl, const hipComplex *d, const hipComplex *du, const hipComplex *dw, const hipComplex *x, int batchCount, size_t *pBufferSizeInBytes)#
hipsparseStatus_t hipsparseZgpsvInterleavedBatch_bufferSizeExt(hipsparseHandle_t handle, int algo, int m, const hipDoubleComplex *ds, const hipDoubleComplex *dl, const hipDoubleComplex *d, const hipDoubleComplex *du, const hipDoubleComplex *dw, const hipDoubleComplex *x, int batchCount, size_t *pBufferSizeInBytes)#

Interleaved Batch pentadiagonal solver.

hipsparseXgpsvInterleavedBatch_bufferSizeExt returns the size of the temporary storage buffer in bytes that is required by hipsparseXgpsvInterleavedBatch(). The temporary storage buffer must be allocated by the user.

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

  • algo[in] algorithm to solve the linear system.

  • m[in] size of the pentadiagonal linear system.

  • ds[in] lower diagonal (distance 2) of pentadiagonal system. First two entries must be zero.

  • dl[in] lower diagonal of pentadiagonal system. First entry must be zero.

  • d[in] main diagonal of pentadiagonal system.

  • du[in] upper diagonal of pentadiagonal system. Last entry must be zero.

  • dw[in] upper diagonal (distance 2) of pentadiagonal system. Last two entries must be zero.

  • x[in] Dense array of right-hand-sides with dimension batchCount by m.

  • batchCount[in] The number of systems to solve.

  • pBufferSizeInBytes[out] Number of bytes of the temporary storage buffer required.

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, alg, batchCount, ds, dl, d, du, dw, x or pBufferSizeInBytes pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.

hipsparseXgpsvInterleavedBatch()#

hipsparseStatus_t hipsparseSgpsvInterleavedBatch(hipsparseHandle_t handle, int algo, int m, float *ds, float *dl, float *d, float *du, float *dw, float *x, int batchCount, void *pBuffer)#
hipsparseStatus_t hipsparseDgpsvInterleavedBatch(hipsparseHandle_t handle, int algo, int m, double *ds, double *dl, double *d, double *du, double *dw, double *x, int batchCount, void *pBuffer)#
hipsparseStatus_t hipsparseCgpsvInterleavedBatch(hipsparseHandle_t handle, int algo, int m, hipComplex *ds, hipComplex *dl, hipComplex *d, hipComplex *du, hipComplex *dw, hipComplex *x, int batchCount, void *pBuffer)#
hipsparseStatus_t hipsparseZgpsvInterleavedBatch(hipsparseHandle_t handle, int algo, int m, hipDoubleComplex *ds, hipDoubleComplex *dl, hipDoubleComplex *d, hipDoubleComplex *du, hipDoubleComplex *dw, hipDoubleComplex *x, int batchCount, void *pBuffer)#

Interleaved Batch pentadiagonal solver.

hipsparseXgpsvInterleavedBatch solves a batched pentadiagonal linear system

Note

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

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

  • algo[in] algorithm to solve the linear system.

  • m[in] size of the pentadiagonal linear system.

  • ds[inout] lower diagonal (distance 2) of pentadiagonal system. First two entries must be zero.

  • dl[inout] lower diagonal of pentadiagonal system. First entry must be zero.

  • d[inout] main diagonal of pentadiagonal system.

  • du[inout] upper diagonal of pentadiagonal system. Last entry must be zero.

  • dw[inout] upper diagonal (distance 2) of pentadiagonal system. Last two entries must be zero.

  • x[inout] Dense array of right-hand-sides with dimension batchCount by m.

  • batchCount[in] The number of systems to solve.

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

Return values:
  • HIPSPARSE_STATUS_SUCCESS – the operation completed successfully.

  • HIPSPARSE_STATUS_INVALID_VALUEhandle, m, alg, batchCount, ds, dl, d, du, dw, x or pBuffer pointer is invalid.

  • HIPSPARSE_STATUS_INTERNAL_ERROR – an internal error occurred.