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.
rocsparse_bsric0_zero_pivot()#
- 
rocsparse_status rocsparse_bsric0_zero_pivot(rocsparse_handle handle, rocsparse_mat_info info, rocsparse_int *position)#
- Incomplete Cholesky factorization with 0 fill-ins and no pivoting using BSR storage format. - rocsparse_bsric0_zero_pivotreturns rocsparse_status_zero_pivot, if either a structural or numerical zero has been found during rocsparse_sbsric0(), rocsparse_dbsric0(), rocsparse_cbsric0() or rocsparse_zbsric0() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored in- position, using same index base as the BSR matrix.- positioncan be in host or device memory. If no zero pivot has been found,- positionis set to -1 and rocsparse_status_success is returned instead.- Note - If a zero pivot is found, - position=jmeans 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 - rocsparse_bsric0_zero_pivotis a blocking function. It might influence performance negatively.- Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_pointer – - infoor- positionpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_zero_pivot – zero pivot has been found. 
 
 
rocsparse_bsric0_buffer_size()#
- 
rocsparse_status rocsparse_sbsric0_buffer_size(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, size_t *buffer_size)#
- 
rocsparse_status rocsparse_dbsric0_buffer_size(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, size_t *buffer_size)#
- 
rocsparse_status rocsparse_cbsric0_buffer_size(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, size_t *buffer_size)#
- 
rocsparse_status rocsparse_zbsric0_buffer_size(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, size_t *buffer_size)#
- Incomplete Cholesky factorization with 0 fill-ins and no pivoting using BSR storage format. - rocsparse_bsric0_buffer_sizereturns the size of the temporary storage buffer that is required by rocsparse_sbsric0_analysis(), rocsparse_dbsric0_analysis(), rocsparse_cbsric0_analysis(), rocsparse_zbsric0_analysis(), rocsparse_sbsric0(), rocsparse_dbsric0(), rocsparse_sbsric0() and rocsparse_dbsric0(). The temporary storage buffer must be allocated by the user. The size of the temporary storage buffer is identical to the size returned by rocsparse_sbsrsv_buffer_size(), rocsparse_dbsrsv_buffer_size(), rocsparse_cbsrsv_buffer_size(), rocsparse_zbsrsv_buffer_size(), rocsparse_sbsrilu0_buffer_size(), rocsparse_dbsrilu0_buffer_size(), rocsparse_cbsrilu0_buffer_size() and rocsparse_zbsrilu0_buffer_size() if the matrix sparsity pattern is identical. The user allocated buffer can thus be shared between subsequent calls to those functions.- Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- dir – [in] direction that specifies whether to count nonzero elements by rocsparse_direction_row or by rocsparse_direction_row. 
- mb – [in] number of block rows in the sparse BSR matrix. 
- nnzb – [in] number of non-zero block entries of the sparse BSR matrix. 
- descr – [in] descriptor of the sparse BSR matrix. 
- bsr_val – [in] array of length - nnzb*block_dim*block_dimcontaining the values of the sparse BSR matrix.
- bsr_row_ptr – [in] array of - mb+1elements that point to the start of every block row of the sparse BSR matrix.
- bsr_col_ind – [in] array of - nnzbelements containing the block column indices of the sparse BSR matrix.
- block_dim – [in] the block dimension of the BSR matrix. Between 1 and m where - m=mb*block_dim.
- info – [out] structure that holds the information collected during the analysis step. 
- buffer_size – [out] number of bytes of the temporary storage buffer required by rocsparse_sbsric0_analysis(), rocsparse_dbsric0_analysis(), rocsparse_cbsric0_analysis(), rocsparse_zbsric0_analysis(), rocsparse_sbsric0(), rocsparse_dbsric0(), rocsparse_cbsric0() and rocsparse_zbsric0(). 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mb,- nnzb, or- block_dimis invalid.
- rocsparse_status_invalid_pointer – - descr,- bsr_val,- bsr_row_ptr,- bsr_col_ind,- infoor- buffer_sizepointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_not_implemented – rocsparse_matrix_type != rocsparse_matrix_type_general. 
 
 
rocsparse_bsric0_analysis()#
- 
rocsparse_status rocsparse_sbsric0_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- 
rocsparse_status rocsparse_dbsric0_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- 
rocsparse_status rocsparse_cbsric0_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- 
rocsparse_status rocsparse_zbsric0_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- Incomplete Cholesky factorization with 0 fill-ins and no pivoting using BSR storage format. - rocsparse_bsric0_analysisperforms the analysis step for rocsparse_sbsric0() rocsparse_dbsric0(), rocsparse_cbsric0(), and rocsparse_zbsric0(). It is expected that this function will be executed only once for a given matrix and particular operation type. The analysis meta data can be cleared by rocsparse_bsric0_clear().- rocsparse_bsric0_analysiscan share its meta data with rocsparse_sbsrilu0_analysis(), rocsparse_dbsrilu0_analysis(), rocsparse_cbsrilu0_analysis(), rocsparse_zbsrilu0_analysis(), rocsparse_sbsrsv_analysis(), rocsparse_dbsrsv_analysis(), rocsparse_cbsrsv_analysis(), rocsparse_zbsrsv_analysis(), rocsparse_sbsrsm_analysis(), rocsparse_dbsrsm_analysis(), rocsparse_cbsrsm_analysis() and rocsparse_zbsrsm_analysis(). Selecting rocsparse_analysis_policy_reuse policy can greatly improve computation performance of meta data. However, the user need to make sure that the sparsity pattern remains unchanged. If this cannot be assured, rocsparse_analysis_policy_force has to be used.- Note - If the matrix sparsity pattern changes, the gathered information will become invalid. - Note - This function is blocking with respect to the host. - Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- dir – [in] direction that specified whether to count nonzero elements by rocsparse_direction_row or by rocsparse_direction_row. 
- mb – [in] number of block rows in the sparse BSR matrix. 
- nnzb – [in] number of non-zero block entries of the sparse BSR matrix. 
- descr – [in] descriptor of the sparse BSR matrix. 
- bsr_val – [in] array of length - nnzb*block_dim*block_dimcontaining the values of the sparse BSR matrix.
- bsr_row_ptr – [in] array of - mb+1elements that point to the start of every block row of the sparse BSR matrix.
- bsr_col_ind – [in] array of - nnzbelements containing the block column indices of the sparse BSR matrix.
- block_dim – [in] the block dimension of the BSR matrix. Between 1 and m where - m=mb*block_dim.
- info – [out] structure that holds the information collected during the analysis step. 
- analysis – [in] rocsparse_analysis_policy_reuse or rocsparse_analysis_policy_force. 
- solve – [in] rocsparse_solve_policy_auto. 
- temp_buffer – [in] temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mb,- nnzb, or- block_dimis invalid.
- rocsparse_status_invalid_pointer – - descr,- bsr_val,- bsr_row_ptr,- bsr_col_ind,- infoor- temp_bufferpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_not_implemented – rocsparse_matrix_type != rocsparse_matrix_type_general. 
 
 
rocsparse_bsric0()#
- 
rocsparse_status rocsparse_sbsric0(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- 
rocsparse_status rocsparse_dbsric0(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- 
rocsparse_status rocsparse_cbsric0(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- 
rocsparse_status rocsparse_zbsric0(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- Incomplete Cholesky factorization with 0 fill-ins and no pivoting using BSR storage format. - rocsparse_bsric0computes 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 \]- rocsparse_bsric0requires a user allocated temporary buffer. Its size is returned by rocsparse_sbsric0_buffer_size(), rocsparse_dbsric0_buffer_size(), rocsparse_cbsric0_buffer_size() or rocsparse_zbsric0_buffer_size(). Furthermore, analysis meta data is required. It can be obtained by rocsparse_sbsric0_analysis(), rocsparse_dbsric0_analysis(), rocsparse_cbsric0_analysis() or rocsparse_zbsric0_analysis().- rocsparse_bsric0reports the first zero pivot (either numerical or structural zero). The zero pivot status can be obtained by calling rocsparse_bsric0_zero_pivot().- Example
- Consider the sparse \(m \times m\) matrix \(A\), stored in BSR storage format. The following example computes the incomplete Cholesky factorization \(M \approx LL^T\) and solves the preconditioned system \(My = x\). - // Create rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); // Create matrix descriptor for M rocsparse_mat_descr descr_M; rocsparse_create_mat_descr(&descr_M); // Create matrix descriptor for L rocsparse_mat_descr descr_L; rocsparse_create_mat_descr(&descr_L); rocsparse_set_mat_fill_mode(descr_L, rocsparse_fill_mode_lower); rocsparse_set_mat_diag_type(descr_L, rocsparse_diag_type_unit); // Create matrix descriptor for L' rocsparse_mat_descr descr_Lt; rocsparse_create_mat_descr(&descr_Lt); rocsparse_set_mat_fill_mode(descr_Lt, rocsparse_fill_mode_upper); rocsparse_set_mat_diag_type(descr_Lt, rocsparse_diag_type_non_unit); // Create matrix info structure rocsparse_mat_info info; rocsparse_create_mat_info(&info); // Obtain required buffer size size_t buffer_size_M; size_t buffer_size_L; size_t buffer_size_Lt; rocsparse_dbsric0_buffer_size(handle, rocsparse_direction_row, mb, nnzb, descr_M, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, &buffer_size_M); rocsparse_dbsrsv_buffer_size(handle, rocsparse_direction_row, rocsparse_operation_none, mb, nnzb, descr_L, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, &buffer_size_L); rocsparse_dbsrsv_buffer_size(handle, rocsparse_direction_row, rocsparse_operation_transpose, mb, nnzb, descr_Lt, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, &buffer_size_Lt); size_t buffer_size = max(buffer_size_M, max(buffer_size_L, buffer_size_Lt)); // Allocate temporary buffer void* temp_buffer; hipMalloc(&temp_buffer, buffer_size); // Perform analysis steps, using rocsparse_analysis_policy_reuse to improve // computation performance rocsparse_dbsric0_analysis(handle, rocsparse_direction_row, mb, nnzb, descr_M, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); rocsparse_dbsrsv_analysis(handle, rocsparse_direction_row, rocsparse_operation_none, mb, nnzb, descr_L, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); rocsparse_dbsrsv_analysis(handle, rocsparse_direction_row, rocsparse_operation_transpose, mb, nnzb, descr_Lt, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); // Check for zero pivot rocsparse_int position; if(rocsparse_status_zero_pivot == rocsparse_bsric0_zero_pivot(handle, info, &position)) { printf("A has structural zero at A(%d,%d)\n", position, position); } // Compute incomplete Cholesky factorization M = LL' rocsparse_dbsric0(handle, rocsparse_direction_row, mb, nnzb, descr_M, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, rocsparse_solve_policy_auto, temp_buffer); // Check for zero pivot if(rocsparse_status_zero_pivot == rocsparse_bsric0_zero_pivot(handle, info, &position)) { printf("L has structural and/or numerical zero at L(%d,%d)\n", position, position); } // Solve Lz = x rocsparse_dbsrsv_solve(handle, rocsparse_direction_row, rocsparse_operation_none, mb, nnzb, &alpha, descr_L, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, x, z, rocsparse_solve_policy_auto, temp_buffer); // Solve L'y = z rocsparse_dbsrsv_solve(handle, rocsparse_direction_row, rocsparse_operation_transpose, mb, nnzb, &alpha, descr_Lt, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, z, y, rocsparse_solve_policy_auto, temp_buffer); // Clean up hipFree(temp_buffer); rocsparse_destroy_mat_info(info); rocsparse_destroy_mat_descr(descr_M); rocsparse_destroy_mat_descr(descr_L); rocsparse_destroy_mat_descr(descr_Lt); rocsparse_destroy_handle(handle); 
 - Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- dir – [in] direction that specified whether to count nonzero elements by rocsparse_direction_row or by rocsparse_direction_row. 
- mb – [in] number of block rows in the sparse BSR matrix. 
- nnzb – [in] number of non-zero block entries of the sparse BSR matrix. 
- descr – [in] descriptor of the sparse BSR matrix. 
- bsr_val – [inout] array of length - nnzb*block_dim*block_dimcontaining the values of the sparse BSR matrix.
- bsr_row_ptr – [in] array of - mb+1elements that point to the start of every block row of the sparse BSR matrix.
- bsr_col_ind – [in] array of - nnzbelements containing the block column indices of the sparse BSR matrix.
- block_dim – [in] the block dimension of the BSR matrix. Between 1 and m where - m=mb*block_dim.
- info – [in] structure that holds the information collected during the analysis step. 
- policy – [in] rocsparse_solve_policy_auto. 
- temp_buffer – [in] temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mb,- nnzb, or- block_dimis invalid.
- rocsparse_status_invalid_pointer – - descr,- bsr_val,- bsr_row_ptror- bsr_col_indpointer is invalid.
- rocsparse_status_arch_mismatch – the device is not supported. 
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_not_implemented – rocsparse_matrix_type != rocsparse_matrix_type_general. 
 
 
rocsparse_bsric0_clear()#
- 
rocsparse_status rocsparse_bsric0_clear(rocsparse_handle handle, rocsparse_mat_info info)#
- Incomplete Cholesky factorization with 0 fill-ins and no pivoting using BSR storage format. - rocsparse_bsric0_cleardeallocates all memory that was allocated by rocsparse_sbsric0_analysis(), rocsparse_dbsric0_analysis(), rocsparse_cbsric0_analysis() or rocsparse_zbsric0_analysis(). This is especially useful, if memory is an issue and the analysis data is not required for further computation.- Note - Calling - rocsparse_bsric0_clearis optional. All allocated resources will be cleared, when the opaque rocsparse_mat_info struct is destroyed using rocsparse_destroy_mat_info().- Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- info – [inout] structure that holds the information collected during the analysis step. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_pointer – - infopointer is invalid.
- rocsparse_status_memory_error – the buffer holding the meta data could not be deallocated. 
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_bsrilu0_zero_pivot()#
- 
rocsparse_status rocsparse_bsrilu0_zero_pivot(rocsparse_handle handle, rocsparse_mat_info info, rocsparse_int *position)#
- Incomplete LU factorization with 0 fill-ins and no pivoting using BSR storage format. - rocsparse_bsrilu0_zero_pivotreturns rocsparse_status_zero_pivot, if either a structural or numerical zero has been found during rocsparse_sbsrilu0(), rocsparse_dbsrilu0(), rocsparse_cbsrilu0() or rocsparse_zbsrilu0() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored in- position, using same index base as the BSR matrix.- positioncan be in host or device memory. If no zero pivot has been found,- positionis set to -1 and rocsparse_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 - rocsparse_bsrilu0_zero_pivotis a blocking function. It might influence performance negatively.- Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_pointer – - infoor- positionpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_zero_pivot – zero pivot has been found. 
 
 
rocsparse_bsrilu0_numeric_boost()#
- 
rocsparse_status rocsparse_sbsrilu0_numeric_boost(rocsparse_handle handle, rocsparse_mat_info info, int enable_boost, const float *boost_tol, const float *boost_val)#
- 
rocsparse_status rocsparse_dbsrilu0_numeric_boost(rocsparse_handle handle, rocsparse_mat_info info, int enable_boost, const double *boost_tol, const double *boost_val)#
- 
rocsparse_status rocsparse_cbsrilu0_numeric_boost(rocsparse_handle handle, rocsparse_mat_info info, int enable_boost, const float *boost_tol, const rocsparse_float_complex *boost_val)#
- 
rocsparse_status rocsparse_zbsrilu0_numeric_boost(rocsparse_handle handle, rocsparse_mat_info info, int enable_boost, const double *boost_tol, const rocsparse_double_complex *boost_val)#
- Incomplete LU factorization with 0 fill-ins and no pivoting using BSR storage format. - rocsparse_bsrilu0_numeric_boostenables the user to replace a numerical value in an incomplete LU factorization.- tolis 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_boostto 1 or disabled by setting- enable_boostto 0.- Note - toland- boost_valcan be in host or device memory.- Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- info – [in] structure that holds the information collected during the analysis step. 
- enable_boost – [in] enable/disable numeric boost. 
- 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:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_pointer – - info,- tolor- boost_valpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_bsrilu0_buffer_size()#
- 
rocsparse_status rocsparse_sbsrilu0_buffer_size(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, size_t *buffer_size)#
- 
rocsparse_status rocsparse_dbsrilu0_buffer_size(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, size_t *buffer_size)#
- 
rocsparse_status rocsparse_cbsrilu0_buffer_size(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, size_t *buffer_size)#
- 
rocsparse_status rocsparse_zbsrilu0_buffer_size(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, size_t *buffer_size)#
- Incomplete LU factorization with 0 fill-ins and no pivoting using BSR storage format. - rocsparse_bsrilu0_buffer_sizereturns the size of the temporary storage buffer that is required by rocsparse_sbsrilu0_analysis(), rocsparse_dbsrilu0_analysis(), rocsparse_cbsrilu0_analysis(), rocsparse_zbsrilu0_analysis(), rocsparse_sbsrilu0(), rocsparse_dbsrilu0(), rocsparse_sbsrilu0() and rocsparse_dbsrilu0(). The temporary storage buffer must be allocated by the user. The size of the temporary storage buffer is identical to the size returned by rocsparse_sbsrsv_buffer_size(), rocsparse_dbsrsv_buffer_size(), rocsparse_cbsrsv_buffer_size(), rocsparse_zbsrsv_buffer_size(), rocsparse_sbsric0_buffer_size(), rocsparse_dbsric0_buffer_size(), rocsparse_cbsric0_buffer_size() and rocsparse_zbsric0_buffer_size() if the matrix sparsity pattern is identical. The user allocated buffer can thus be shared between subsequent calls to those functions.- Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- dir – [in] direction that specifies whether to count nonzero elements by rocsparse_direction_row or by rocsparse_direction_row. 
- mb – [in] number of block rows in the sparse BSR matrix. 
- nnzb – [in] number of non-zero block entries of the sparse BSR matrix. 
- descr – [in] descriptor of the sparse BSR matrix. 
- bsr_val – [in] array of length - nnzb*block_dim*block_dimcontaining the values of the sparse BSR matrix.
- bsr_row_ptr – [in] array of - mb+1elements that point to the start of every block row of the sparse BSR matrix.
- bsr_col_ind – [in] array of - nnzbelements containing the block column indices of the sparse BSR matrix.
- block_dim – [in] the block dimension of the BSR matrix. Between 1 and m where - m=mb*block_dim.
- info – [out] structure that holds the information collected during the analysis step. 
- buffer_size – [out] number of bytes of the temporary storage buffer required by rocsparse_sbsrilu0_analysis(), rocsparse_dbsrilu0_analysis(), rocsparse_cbsrilu0_analysis(), rocsparse_zbsrilu0_analysis(), rocsparse_sbsrilu0(), rocsparse_dbsrilu0(), rocsparse_cbsrilu0() and rocsparse_zbsrilu0(). 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mb,- nnzb, or- block_dimis invalid.
- rocsparse_status_invalid_pointer – - descr,- bsr_val,- bsr_row_ptr,- bsr_col_ind,- infoor- buffer_sizepointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_not_implemented – rocsparse_matrix_type != rocsparse_matrix_type_general. 
 
 
rocsparse_bsrilu0_analysis()#
- 
rocsparse_status rocsparse_sbsrilu0_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- 
rocsparse_status rocsparse_dbsrilu0_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- 
rocsparse_status rocsparse_cbsrilu0_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- 
rocsparse_status rocsparse_zbsrilu0_analysis(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, const rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- Incomplete LU factorization with 0 fill-ins and no pivoting using BSR storage format. - rocsparse_bsrilu0_analysisperforms the analysis step for rocsparse_sbsrilu0() rocsparse_dbsrilu0(), rocsparse_cbsrilu0(), and rocsparse_zbsrilu0(). It is expected that this function will be executed only once for a given matrix. The analysis meta data can be cleared by rocsparse_bsrilu0_clear().- rocsparse_bsrilu0_analysiscan share its meta data with rocsparse_sbsric0_analysis(), rocsparse_dbsric0_analysis(), rocsparse_cbsric0_analysis(), rocsparse_zbsric0_analysis(), rocsparse_sbsrsv_analysis(), rocsparse_dbsrsv_analysis(), rocsparse_cbsrsv_analysis(), rocsparse_zbsrsv_analysis(), rocsparse_sbsrsm_analysis(), rocsparse_dbsrsm_analysis(), rocsparse_cbsrsm_analysis() and rocsparse_zbsrsm_analysis(). Selecting rocsparse_analysis_policy_reuse policy can greatly improve computation performance of meta data. However, the user need to make sure that the sparsity pattern remains unchanged. If this cannot be assured, rocsparse_analysis_policy_force has to be used.- Note - If the matrix sparsity pattern changes, the gathered information will become invalid. - Note - This function is blocking with respect to the host. - Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- dir – [in] direction that specified whether to count nonzero elements by rocsparse_direction_row or by rocsparse_direction_row. 
- mb – [in] number of block rows in the sparse BSR matrix. 
- nnzb – [in] number of non-zero block entries of the sparse BSR matrix. 
- descr – [in] descriptor of the sparse BSR matrix. 
- bsr_val – [in] array of length - nnzb*block_dim*block_dimcontaining the values of the sparse BSR matrix.
- bsr_row_ptr – [in] array of - mb+1elements that point to the start of every block row of the sparse BSR matrix.
- bsr_col_ind – [in] array of - nnzbelements containing the block column indices of the sparse BSR matrix.
- block_dim – [in] the block dimension of the BSR matrix. Between 1 and m where - m=mb*block_dim.
- info – [out] structure that holds the information collected during the analysis step. 
- analysis – [in] rocsparse_analysis_policy_reuse or rocsparse_analysis_policy_force. 
- solve – [in] rocsparse_solve_policy_auto. 
- temp_buffer – [in] temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mb,- nnzb, or- block_dimis invalid.
- rocsparse_status_invalid_pointer – - descr,- bsr_val,- bsr_row_ptr,- bsr_col_ind,- infoor- temp_bufferpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_not_implemented – rocsparse_matrix_type != rocsparse_matrix_type_general. 
 
 
rocsparse_bsrilu0()#
- 
rocsparse_status rocsparse_sbsrilu0(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, float *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- 
rocsparse_status rocsparse_dbsrilu0(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, double *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- 
rocsparse_status rocsparse_cbsrilu0(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, rocsparse_float_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- 
rocsparse_status rocsparse_zbsrilu0(rocsparse_handle handle, rocsparse_direction dir, rocsparse_int mb, rocsparse_int nnzb, const rocsparse_mat_descr descr, rocsparse_double_complex *bsr_val, const rocsparse_int *bsr_row_ptr, const rocsparse_int *bsr_col_ind, rocsparse_int block_dim, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- Incomplete LU factorization with 0 fill-ins and no pivoting using BSR storage format. - rocsparse_bsrilu0computes 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 \]- rocsparse_bsrilu0requires a user allocated temporary buffer. Its size is returned by rocsparse_sbsrilu0_buffer_size(), rocsparse_dbsrilu0_buffer_size(), rocsparse_cbsrilu0_buffer_size() or rocsparse_zbsrilu0_buffer_size(). Furthermore, analysis meta data is required. It can be obtained by rocsparse_sbsrilu0_analysis(), rocsparse_dbsrilu0_analysis(), rocsparse_cbsrilu0_analysis() or rocsparse_zbsrilu0_analysis().- rocsparse_bsrilu0reports the first zero pivot (either numerical or structural zero). The zero pivot status can be obtained by calling rocsparse_bsrilu0_zero_pivot().- Example
- Consider the sparse \(m \times m\) matrix \(A\), stored in BSR storage format. The following example computes the incomplete LU factorization \(M \approx LU\) and solves the preconditioned system \(My = x\). - // Create rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); // Create matrix descriptor for M rocsparse_mat_descr descr_M; rocsparse_create_mat_descr(&descr_M); // Create matrix descriptor for L rocsparse_mat_descr descr_L; rocsparse_create_mat_descr(&descr_L); rocsparse_set_mat_fill_mode(descr_L, rocsparse_fill_mode_lower); rocsparse_set_mat_diag_type(descr_L, rocsparse_diag_type_unit); // Create matrix descriptor for U rocsparse_mat_descr descr_U; rocsparse_create_mat_descr(&descr_U); rocsparse_set_mat_fill_mode(descr_U, rocsparse_fill_mode_upper); rocsparse_set_mat_diag_type(descr_U, rocsparse_diag_type_non_unit); // Create matrix info structure rocsparse_mat_info info; rocsparse_create_mat_info(&info); // Obtain required buffer size size_t buffer_size_M; size_t buffer_size_L; size_t buffer_size_U; rocsparse_dbsrilu0_buffer_size(handle, rocsparse_direction_row, mb, nnzb, descr_M, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, &buffer_size_M); rocsparse_dbsrsv_buffer_size(handle, rocsparse_direction_row, rocsparse_operation_none, mb, nnzb, descr_L, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, &buffer_size_L); rocsparse_dbsrsv_buffer_size(handle, rocsparse_direction_row, rocsparse_operation_none, mb, nnzb, descr_U, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, &buffer_size_U); size_t buffer_size = max(buffer_size_M, max(buffer_size_L, buffer_size_U)); // Allocate temporary buffer void* temp_buffer; hipMalloc(&temp_buffer, buffer_size); // Perform analysis steps, using rocsparse_analysis_policy_reuse to improve // computation performance rocsparse_dbsrilu0_analysis(handle, rocsparse_direction_row, mb, nnzb, descr_M, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); rocsparse_dbsrsv_analysis(handle, rocsparse_direction_row, rocsparse_operation_none, mb, nnzb, descr_L, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); rocsparse_dbsrsv_analysis(handle, rocsparse_direction_row, rocsparse_operation_none, mb, nnzb, descr_U, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); // Check for zero pivot rocsparse_int position; if(rocsparse_status_zero_pivot == rocsparse_bsrilu0_zero_pivot(handle, info, &position)) { printf("A has structural zero at A(%d,%d)\n", position, position); } // Compute incomplete LU factorization M = LU rocsparse_dbsrilu0(handle, rocsparse_direction_row, mb, nnzb, descr_M, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, rocsparse_solve_policy_auto, temp_buffer); // Check for zero pivot if(rocsparse_status_zero_pivot == rocsparse_bsrilu0_zero_pivot(handle, info, &position)) { printf("L has structural and/or numerical zero at L(%d,%d)\n", position, position); } // Solve Lz = x rocsparse_dbsrsv_solve(handle, rocsparse_direction_row, rocsparse_operation_none, mb, nnzb, &alpha, descr_L, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, x, z, rocsparse_solve_policy_auto, temp_buffer); // Solve Uy = z rocsparse_dbsrsv_solve(handle, rocsparse_direction_row, rocsparse_operation_none, mb, nnzb, &alpha, descr_U, bsr_val, bsr_row_ptr, bsr_col_ind, block_dim, info, z, y, rocsparse_solve_policy_auto, temp_buffer); // Clean up hipFree(temp_buffer); rocsparse_destroy_mat_info(info); rocsparse_destroy_mat_descr(descr_M); rocsparse_destroy_mat_descr(descr_L); rocsparse_destroy_mat_descr(descr_U); rocsparse_destroy_handle(handle); 
 - Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- dir – [in] direction that specified whether to count nonzero elements by rocsparse_direction_row or by rocsparse_direction_row. 
- mb – [in] number of block rows in the sparse BSR matrix. 
- nnzb – [in] number of non-zero block entries of the sparse BSR matrix. 
- descr – [in] descriptor of the sparse BSR matrix. 
- bsr_val – [inout] array of length - nnzb*block_dim*block_dimcontaining the values of the sparse BSR matrix.
- bsr_row_ptr – [in] array of - mb+1elements that point to the start of every block row of the sparse BSR matrix.
- bsr_col_ind – [in] array of - nnzbelements containing the block column indices of the sparse BSR matrix.
- block_dim – [in] the block dimension of the BSR matrix. Between 1 and m where - m=mb*block_dim.
- info – [in] structure that holds the information collected during the analysis step. 
- policy – [in] rocsparse_solve_policy_auto. 
- temp_buffer – [in] temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mb,- nnzb, or- block_dimis invalid.
- rocsparse_status_invalid_pointer – - descr,- bsr_val,- bsr_row_ptror- bsr_col_indpointer is invalid.
- rocsparse_status_arch_mismatch – the device is not supported. 
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_not_implemented – rocsparse_matrix_type != rocsparse_matrix_type_general. 
 
 
rocsparse_bsrilu0_clear()#
- 
rocsparse_status rocsparse_bsrilu0_clear(rocsparse_handle handle, rocsparse_mat_info info)#
- Incomplete LU factorization with 0 fill-ins and no pivoting using BSR storage format. - rocsparse_bsrilu0_cleardeallocates all memory that was allocated by rocsparse_sbsrilu0_analysis(), rocsparse_dbsrilu0_analysis(), rocsparse_cbsrilu0_analysis() or rocsparse_zbsrilu0_analysis(). This is especially useful, if memory is an issue and the analysis data is not required for further computation.- Note - Calling - rocsparse_bsrilu0_clearis optional. All allocated resources will be cleared, when the opaque rocsparse_mat_info struct is destroyed using rocsparse_destroy_mat_info().- Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- info – [inout] structure that holds the information collected during the analysis step. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_pointer – - infopointer is invalid.
- rocsparse_status_memory_error – the buffer holding the meta data could not be deallocated. 
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_csric0_zero_pivot()#
- 
rocsparse_status rocsparse_csric0_zero_pivot(rocsparse_handle handle, rocsparse_mat_info info, rocsparse_int *position)#
- Incomplete Cholesky factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csric_zero_pivotreturns rocsparse_status_zero_pivot, if either a structural or numerical zero has been found during rocsparse_scsric0() or rocsparse_dcsric0() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored in- position, using same index base as the CSR matrix.- positioncan be in host or device memory. If no zero pivot has been found,- positionis set to -1 and rocsparse_status_success is returned instead.- Note - rocsparse_csric0_zero_pivotis a blocking function. It might influence performance negatively.- Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_pointer – - infoor- positionpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_zero_pivot – zero pivot has been found. 
 
 
rocsparse_csric0_buffer_size()#
- 
rocsparse_status rocsparse_scsric0_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
- 
rocsparse_status rocsparse_dcsric0_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
- 
rocsparse_status rocsparse_ccsric0_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
- 
rocsparse_status rocsparse_zcsric0_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
- Incomplete Cholesky factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csric0_buffer_sizereturns the size of the temporary storage buffer that is required by rocsparse_scsric0_analysis(), rocsparse_dcsric0_analysis(), rocsparse_scsric0() and rocsparse_dcsric0(). The temporary storage buffer must be allocated by the user. The size of the temporary storage buffer is identical to the size returned by rocsparse_scsrsv_buffer_size(), rocsparse_dcsrsv_buffer_size(), rocsparse_scsrilu0_buffer_size() and rocsparse_dcsrilu0_buffer_size() if the matrix sparsity pattern is identical. The user allocated buffer can thus be shared between subsequent calls to those functions.- Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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. 
- descr – [in] descriptor of the sparse CSR matrix. 
- csr_val – [in] array of - nnzelements of the sparse CSR matrix.
- csr_row_ptr – [in] array of - m+1elements that point to the start of every row of the sparse CSR matrix.
- csr_col_ind – [in] array of - nnzelements containing the column indices of the sparse CSR matrix.
- info – [out] structure that holds the information collected during the analysis step. 
- buffer_size – [out] number of bytes of the temporary storage buffer required by rocsparse_scsric0_analysis(), rocsparse_dcsric0_analysis(), rocsparse_scsric0() and rocsparse_dcsric0(). 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mor- nnzis invalid.
- rocsparse_status_invalid_pointer – - descr,- csr_val,- csr_row_ptr,- csr_col_ind,- infoor- buffer_sizepointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_not_implemented – - trans!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
 
 
rocsparse_csric0_analysis()#
- 
rocsparse_status rocsparse_scsric0_analysis(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- 
rocsparse_status rocsparse_dcsric0_analysis(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- 
rocsparse_status rocsparse_ccsric0_analysis(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- 
rocsparse_status rocsparse_zcsric0_analysis(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- Incomplete Cholesky factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csric0_analysisperforms the analysis step for rocsparse_scsric0() and rocsparse_dcsric0(). It is expected that this function will be executed only once for a given matrix and particular operation type. The analysis meta data can be cleared by rocsparse_csric0_clear().- rocsparse_csric0_analysiscan share its meta data with rocsparse_scsrilu0_analysis(), rocsparse_dcsrilu0_analysis(), rocsparse_ccsrilu0_analysis(), rocsparse_zcsrilu0_analysis(), rocsparse_scsrsv_analysis(), rocsparse_dcsrsv_analysis(), rocsparse_ccsrsv_analysis(), rocsparse_zcsrsv_analysis(), rocsparse_scsrsm_analysis(), rocsparse_dcsrsm_analysis(), rocsparse_scsrsm_analysis() and rocsparse_dcsrsm_analysis(). Selecting rocsparse_analysis_policy_reuse policy can greatly improve computation performance of meta data. However, the user need to make sure that the sparsity pattern remains unchanged. If this cannot be assured, rocsparse_analysis_policy_force has to be used.- Note - If the matrix sparsity pattern changes, the gathered information will become invalid. - Note - This function is blocking with respect to the host. - Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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. 
- descr – [in] descriptor of the sparse CSR matrix. 
- csr_val – [in] array of - nnzelements of the sparse CSR matrix.
- csr_row_ptr – [in] array of - m+1elements that point to the start of every row of the sparse CSR matrix.
- csr_col_ind – [in] array of - nnzelements containing the column indices of the sparse CSR matrix.
- info – [out] structure that holds the information collected during the analysis step. 
- analysis – [in] rocsparse_analysis_policy_reuse or rocsparse_analysis_policy_force. 
- solve – [in] rocsparse_solve_policy_auto. 
- temp_buffer – [in] temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mor- nnzis invalid.
- rocsparse_status_invalid_pointer – - descr,- csr_val,- csr_row_ptr,- csr_col_ind,- infoor- temp_bufferpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_not_implemented – - trans!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
 
 
rocsparse_csric0()#
- 
rocsparse_status rocsparse_scsric0(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- 
rocsparse_status rocsparse_dcsric0(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- 
rocsparse_status rocsparse_ccsric0(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- 
rocsparse_status rocsparse_zcsric0(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- Incomplete Cholesky factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csric0computes 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 \]- rocsparse_csric0requires a user allocated temporary buffer. Its size is returned by rocsparse_scsric0_buffer_size() or rocsparse_dcsric0_buffer_size(). Furthermore, analysis meta data is required. It can be obtained by rocsparse_scsric0_analysis() or rocsparse_dcsric0_analysis().- rocsparse_csric0reports the first zero pivot (either numerical or structural zero). The zero pivot status can be obtained by calling rocsparse_csric0_zero_pivot().- Example
- Consider the sparse \(m \times m\) matrix \(A\), stored in CSR storage format. The following example computes the incomplete Cholesky factorization \(M \approx LL^T\) and solves the preconditioned system \(My = x\). - // Create rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); // Create matrix descriptor for M rocsparse_mat_descr descr_M; rocsparse_create_mat_descr(&descr_M); // Create matrix descriptor for L rocsparse_mat_descr descr_L; rocsparse_create_mat_descr(&descr_L); rocsparse_set_mat_fill_mode(descr_L, rocsparse_fill_mode_lower); rocsparse_set_mat_diag_type(descr_L, rocsparse_diag_type_unit); // Create matrix descriptor for L' rocsparse_mat_descr descr_Lt; rocsparse_create_mat_descr(&descr_Lt); rocsparse_set_mat_fill_mode(descr_Lt, rocsparse_fill_mode_upper); rocsparse_set_mat_diag_type(descr_Lt, rocsparse_diag_type_non_unit); // Create matrix info structure rocsparse_mat_info info; rocsparse_create_mat_info(&info); // Obtain required buffer size size_t buffer_size_M; size_t buffer_size_L; size_t buffer_size_Lt; rocsparse_dcsric0_buffer_size(handle, m, nnz, descr_M, csr_val, csr_row_ptr, csr_col_ind, info, &buffer_size_M); rocsparse_dcsrsv_buffer_size(handle, rocsparse_operation_none, m, nnz, descr_L, csr_val, csr_row_ptr, csr_col_ind, info, &buffer_size_L); rocsparse_dcsrsv_buffer_size(handle, rocsparse_operation_transpose, m, nnz, descr_Lt, csr_val, csr_row_ptr, csr_col_ind, info, &buffer_size_Lt); size_t buffer_size = max(buffer_size_M, max(buffer_size_L, buffer_size_Lt)); // Allocate temporary buffer void* temp_buffer; hipMalloc(&temp_buffer, buffer_size); // Perform analysis steps, using rocsparse_analysis_policy_reuse to improve // computation performance rocsparse_dcsric0_analysis(handle, m, nnz, descr_M, csr_val, csr_row_ptr, csr_col_ind, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); rocsparse_dcsrsv_analysis(handle, rocsparse_operation_none, m, nnz, descr_L, csr_val, csr_row_ptr, csr_col_ind, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); rocsparse_dcsrsv_analysis(handle, rocsparse_operation_transpose, m, nnz, descr_Lt, csr_val, csr_row_ptr, csr_col_ind, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); // Check for zero pivot rocsparse_int position; if(rocsparse_status_zero_pivot == rocsparse_csric0_zero_pivot(handle, info, &position)) { printf("A has structural zero at A(%d,%d)\n", position, position); } // Compute incomplete Cholesky factorization M = LL' rocsparse_dcsric0(handle, m, nnz, descr_M, csr_val, csr_row_ptr, csr_col_ind, info, rocsparse_solve_policy_auto, temp_buffer); // Check for zero pivot if(rocsparse_status_zero_pivot == rocsparse_csric0_zero_pivot(handle, info, &position)) { printf("L has structural and/or numerical zero at L(%d,%d)\n", position, position); } // Solve Lz = x rocsparse_dcsrsv_solve(handle, rocsparse_operation_none, m, nnz, &alpha, descr_L, csr_val, csr_row_ptr, csr_col_ind, info, x, z, rocsparse_solve_policy_auto, temp_buffer); // Solve L'y = z rocsparse_dcsrsv_solve(handle, rocsparse_operation_transpose, m, nnz, &alpha, descr_Lt, csr_val, csr_row_ptr, csr_col_ind, info, z, y, rocsparse_solve_policy_auto, temp_buffer); // Clean up hipFree(temp_buffer); rocsparse_destroy_mat_info(info); rocsparse_destroy_mat_descr(descr_M); rocsparse_destroy_mat_descr(descr_L); rocsparse_destroy_mat_descr(descr_Lt); rocsparse_destroy_handle(handle); 
 - Note - The sparse CSR matrix has to be sorted. This can be achieved by calling rocsparse_csrsort(). - Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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. 
- descr – [in] descriptor of the sparse CSR matrix. 
- csr_val – [inout] array of - nnzelements of the sparse CSR matrix.
- csr_row_ptr – [in] array of - m+1elements that point to the start of every row of the sparse CSR matrix.
- csr_col_ind – [in] array of - nnzelements containing the column indices of the sparse CSR matrix.
- info – [in] structure that holds the information collected during the analysis step. 
- policy – [in] rocsparse_solve_policy_auto. 
- temp_buffer – [in] temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mor- nnzis invalid.
- rocsparse_status_invalid_pointer – - descr,- csr_val,- csr_row_ptror- csr_col_indpointer is invalid.
- rocsparse_status_arch_mismatch – the device is not supported. 
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_not_implemented – - trans!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
 
 
rocsparse_csric0_clear()#
- 
rocsparse_status rocsparse_csric0_clear(rocsparse_handle handle, rocsparse_mat_info info)#
- Incomplete Cholesky factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csric0_cleardeallocates all memory that was allocated by rocsparse_scsric0_analysis() or rocsparse_dcsric0_analysis(). This is especially useful, if memory is an issue and the analysis data is not required for further computation.- Note - Calling - rocsparse_csric0_clearis optional. All allocated resources will be cleared, when the opaque rocsparse_mat_info struct is destroyed using rocsparse_destroy_mat_info().- Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- info – [inout] structure that holds the information collected during the analysis step. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_pointer – - infopointer is invalid.
- rocsparse_status_memory_error – the buffer holding the meta data could not be deallocated. 
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_csritilu0_buffer_size()#
- 
rocsparse_status rocsparse_csritilu0_buffer_size(rocsparse_handle handle, rocsparse_itilu0_alg alg, rocsparse_int option, rocsparse_int nmaxiter, rocsparse_int m, rocsparse_int nnz, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_index_base idx_base, rocsparse_datatype datatype, size_t *buffer_size)#
- Iterative Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csritilu0_buffer_sizecomputes the size in bytes of the buffer that has to be allocated by the user.- Note - The sparse CSR matrix has to be sorted. This can be achieved by calling rocsparse_csrsort(). - Note - This function is blocking with respect to the host. - Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- alg – [in] algorithm to use, rocsparse_itilu0_alg 
- option – [in] combination of enumeration values from rocsparse_itilu0_option. 
- nmaxiter – [in] maximum number of iterations. 
- m – [in] number of rows of the sparse CSR matrix. 
- nnz – [in] number of non-zero entries of the sparse CSR matrix. 
- csr_row_ptr – [in] array of - m+1elements that point to the start of every row of the sparse CSR matrix.
- csr_col_ind – [in] array of - nnzelements containing the column indices of the sparse CSR matrix.
- idx_base – [in] rocsparse_index_base_zero or rocsparse_index_base_one. 
- datatype – [in] Type of numerical values, rocsparse_datatype. 
- buffer_size – [out] size of the temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mor- nnzis invalid.
- rocsparse_status_invalid_value – - alg,- baseor datatype is invalid.
- rocsparse_status_invalid_pointer – - csr_row_ptror- csr_col_indpointer is invalid.
- rocsparse_status_zero_pivot – if nnz is zero. 
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_csritilu0_preprocess()#
- 
rocsparse_status rocsparse_csritilu0_preprocess(rocsparse_handle handle, rocsparse_itilu0_alg alg, rocsparse_int option, rocsparse_int nmaxiter, rocsparse_int m, rocsparse_int nnz, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_index_base idx_base, rocsparse_datatype datatype, size_t buffer_size, void *buffer)#
- Iterative Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csritilu0_preprocesscomputes the information required to run rocsparse_scsritilu0_compute, rocsparse_dcsritilu0_compute, rocsparse_ccsritilu0_compute, or rocsparse_zcsritilu0_compute, and stores it in the buffer.- Note - The sparse CSR matrix has to be sorted. This can be achieved by calling rocsparse_csrsort(). - Note - This function is blocking with respect to the host. - Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- alg – [in] algorithm to use, rocsparse_itilu0_alg 
- option – [in] combination of enumeration values from rocsparse_itilu0_option. 
- nmaxiter – [in] maximum number of iterations. 
- m – [in] number of rows of the sparse CSR matrix. 
- nnz – [in] number of non-zero entries of the sparse CSR matrix. 
- csr_row_ptr – [in] array of - m+1elements that point to the start of every row of the sparse CSR matrix.
- csr_col_ind – [in] array of - nnzelements containing the column indices of the sparse CSR matrix.
- idx_base – [in] rocsparse_index_base_zero or rocsparse_index_base_one. 
- datatype – [in] Type of numerical values, rocsparse_datatype. 
- buffer_size – [in] size of the storage buffer allocated by the user. 
- buffer – [in] storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_value – - alg,- baseor datatype is invalid.
- rocsparse_status_invalid_size – - mor- nnzis invalid.
- rocsparse_status_invalid_pointer – - csr_row_ptror- csr_col_indpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_zero_pivot – if missing diagonal element is detected. 
 
 
rocsparse_csritilu0_history()#
- 
rocsparse_status rocsparse_scsritilu0_history(rocsparse_handle handle, rocsparse_itilu0_alg alg, rocsparse_int *niter, float *data, size_t buffer_size, void *buffer)#
- 
rocsparse_status rocsparse_dcsritilu0_history(rocsparse_handle handle, rocsparse_itilu0_alg alg, rocsparse_int *niter, double *data, size_t buffer_size, void *buffer)#
- 
rocsparse_status rocsparse_ccsritilu0_history(rocsparse_handle handle, rocsparse_itilu0_alg alg, rocsparse_int *niter, float *data, size_t buffer_size, void *buffer)#
- 
rocsparse_status rocsparse_zcsritilu0_history(rocsparse_handle handle, rocsparse_itilu0_alg alg, rocsparse_int *niter, double *data, size_t buffer_size, void *buffer)#
- Iterative Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csritilu0_historyfetches convergence history data.- Note - The sparse CSR matrix has to be sorted. This can be achieved by calling rocsparse_csrsort(). - Note - This function is blocking with respect to the host. - Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- alg – [in] algorithm to use, rocsparse_itilu0_alg 
- niter – [out] number of performed iterations. 
- data – [out] norms. 
- buffer_size – [in] size of the buffer allocated by the user. 
- buffer – [in] buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_pointer – - niteror- datais invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_csritilu0_compute()#
- 
rocsparse_status rocsparse_scsritilu0_compute(rocsparse_handle handle, rocsparse_itilu0_alg alg, rocsparse_int option, rocsparse_int *nmaxiter, float tol, rocsparse_int m, rocsparse_int nnz, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, const float *csr_val, float *ilu0, rocsparse_index_base idx_base, size_t buffer_size, void *buffer)#
- 
rocsparse_status rocsparse_dcsritilu0_compute(rocsparse_handle handle, rocsparse_itilu0_alg alg, rocsparse_int option, rocsparse_int *nmaxiter, double tol, rocsparse_int m, rocsparse_int nnz, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, const double *csr_val, double *ilu0, rocsparse_index_base idx_base, size_t buffer_size, void *buffer)#
- 
rocsparse_status rocsparse_ccsritilu0_compute(rocsparse_handle handle, rocsparse_itilu0_alg alg, rocsparse_int option, rocsparse_int *nmaxiter, float tol, rocsparse_int m, rocsparse_int nnz, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, const rocsparse_float_complex *csr_val, rocsparse_float_complex *ilu0, rocsparse_index_base idx_base, size_t buffer_size, void *buffer)#
- 
rocsparse_status rocsparse_zcsritilu0_compute(rocsparse_handle handle, rocsparse_itilu0_alg alg, rocsparse_int option, rocsparse_int *nmaxiter, double tol, rocsparse_int m, rocsparse_int nnz, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, const rocsparse_double_complex *csr_val, rocsparse_double_complex *ilu0, rocsparse_index_base idx_base, size_t buffer_size, void *buffer)#
- Iterative Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csritilu0_computecomputes iteratively 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 \]- rocsparse_csritilu0requires a user allocated temporary buffer. Its size is returned by rocsparse_csritilu0_buffer_size(). Furthermore, analysis meta data is required. It can be obtained by rocsparse_csritlu0_preprocess().- Note - The sparse CSR matrix has to be sorted. This can be achieved by calling rocsparse_csrsort(). - Note - This function is blocking with respect to the host. - Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- alg – [in] algorithm to use, rocsparse_itilu0_alg 
- option – [in] combination of enumeration values from rocsparse_itilu0_option. 
- nmaxiter – [inout] maximum number of iterations. 
- tol – [in] tolerance to use for stopping criteria. 
- m – [in] number of rows of the sparse CSR matrix. 
- nnz – [in] number of non-zero entries of the sparse CSR matrix. 
- csr_row_ptr – [in] array of - m+1elements that point to the start of every row of the sparse CSR matrix.
- csr_col_ind – [in] array of - nnzelements containing the column indices of the sparse CSR matrix.
- csr_val – [inout] array of - nnzelements of the sparse CSR matrix.
- ilu0 – [out] incomplete factorization. 
- idx_base – [in] rocsparse_index_base_zero or rocsparse_index_base_one. 
- buffer_size – [in] size of the storage buffer allocated by the user. 
- buffer – [in] storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_value – - algor- baseis invalid.
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mor- nnzis invalid.
- rocsparse_status_invalid_pointer – - csr_row_ptror- csr_col_indpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_csrilu0_zero_pivot()#
- 
rocsparse_status rocsparse_csrilu0_zero_pivot(rocsparse_handle handle, rocsparse_mat_info info, rocsparse_int *position)#
- Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csrilu0_zero_pivotreturns rocsparse_status_zero_pivot, if either a structural or numerical zero has been found during rocsparse_scsrilu0(), rocsparse_dcsrilu0(), rocsparse_ccsrilu0() or rocsparse_zcsrilu0() computation. The first zero pivot \(j\) at \(A_{j,j}\) is stored in- position, using same index base as the CSR matrix.- positioncan be in host or device memory. If no zero pivot has been found,- positionis set to -1 and rocsparse_status_success is returned instead.- Note - rocsparse_csrilu0_zero_pivotis a blocking function. It might influence performance negatively.- Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_pointer – - infoor- positionpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_zero_pivot – zero pivot has been found. 
 
 
rocsparse_csrilu0_numeric_boost()#
- 
rocsparse_status rocsparse_scsrilu0_numeric_boost(rocsparse_handle handle, rocsparse_mat_info info, int enable_boost, const float *boost_tol, const float *boost_val)#
- 
rocsparse_status rocsparse_dcsrilu0_numeric_boost(rocsparse_handle handle, rocsparse_mat_info info, int enable_boost, const double *boost_tol, const double *boost_val)#
- 
rocsparse_status rocsparse_ccsrilu0_numeric_boost(rocsparse_handle handle, rocsparse_mat_info info, int enable_boost, const float *boost_tol, const rocsparse_float_complex *boost_val)#
- 
rocsparse_status rocsparse_zcsrilu0_numeric_boost(rocsparse_handle handle, rocsparse_mat_info info, int enable_boost, const double *boost_tol, const rocsparse_double_complex *boost_val)#
- Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csrilu0_numeric_boostenables the user to replace a numerical value in an incomplete LU factorization.- tolis 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_boostto 1 or disabled by setting- enable_boostto 0.- Note - toland- boost_valcan be in host or device memory.- Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- info – [in] structure that holds the information collected during the analysis step. 
- enable_boost – [in] enable/disable numeric boost. 
- 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:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_pointer – - info,- tolor- boost_valpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_csrilu0_buffer_size()#
- 
rocsparse_status rocsparse_scsrilu0_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
- 
rocsparse_status rocsparse_dcsrilu0_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
- 
rocsparse_status rocsparse_ccsrilu0_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
- 
rocsparse_status rocsparse_zcsrilu0_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, size_t *buffer_size)#
- Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csrilu0_buffer_sizereturns the size of the temporary storage buffer that is required by rocsparse_scsrilu0_analysis(), rocsparse_dcsrilu0_analysis(), rocsparse_ccsrilu0_analysis(), rocsparse_zcsrilu0_analysis(), rocsparse_scsrilu0(), rocsparse_dcsrilu0(), rocsparse_ccsrilu0() and rocsparse_zcsrilu0(). The temporary storage buffer must be allocated by the user. The size of the temporary storage buffer is identical to the size returned by rocsparse_scsrsv_buffer_size(), rocsparse_dcsrsv_buffer_size(), rocsparse_ccsrsv_buffer_size() and rocsparse_zcsrsv_buffer_size() if the matrix sparsity pattern is identical. The user allocated buffer can thus be shared between subsequent calls to those functions.- Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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. 
- descr – [in] descriptor of the sparse CSR matrix. 
- csr_val – [in] array of - nnzelements of the sparse CSR matrix.
- csr_row_ptr – [in] array of - m+1elements that point to the start of every row of the sparse CSR matrix.
- csr_col_ind – [in] array of - nnzelements containing the column indices of the sparse CSR matrix.
- info – [out] structure that holds the information collected during the analysis step. 
- buffer_size – [out] number of bytes of the temporary storage buffer required by rocsparse_scsrilu0_analysis(), rocsparse_dcsrilu0_analysis(), rocsparse_ccsrilu0_analysis(), rocsparse_zcsrilu0_analysis(), rocsparse_scsrilu0(), rocsparse_dcsrilu0(), rocsparse_ccsrilu0() and rocsparse_zcsrilu0(). 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mor- nnzis invalid.
- rocsparse_status_invalid_pointer – - descr,- csr_val,- csr_row_ptr,- csr_col_ind,- infoor- buffer_sizepointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_not_implemented – - trans!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
 
 
rocsparse_csrilu0_analysis()#
- 
rocsparse_status rocsparse_scsrilu0_analysis(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- 
rocsparse_status rocsparse_dcsrilu0_analysis(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- 
rocsparse_status rocsparse_ccsrilu0_analysis(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- 
rocsparse_status rocsparse_zcsrilu0_analysis(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, const rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_analysis_policy analysis, rocsparse_solve_policy solve, void *temp_buffer)#
- Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csrilu0_analysisperforms the analysis step for rocsparse_scsrilu0(), rocsparse_dcsrilu0(), rocsparse_ccsrilu0() and rocsparse_zcsrilu0(). It is expected that this function will be executed only once for a given matrix and particular operation type. The analysis meta data can be cleared by rocsparse_csrilu0_clear().- rocsparse_csrilu0_analysiscan share its meta data with rocsparse_scsric0_analysis(), rocsparse_dcsric0_analysis(), rocsparse_ccsric0_analysis(), rocsparse_zcsric0_analysis(), rocsparse_scsrsv_analysis(), rocsparse_dcsrsv_analysis(), rocsparse_ccsrsv_analysis(), rocsparse_zcsrsv_analysis(), rocsparse_scsrsm_analysis(), rocsparse_dcsrsm_analysis(), rocsparse_scsrsm_analysis() and rocsparse_dcsrsm_analysis(). Selecting rocsparse_analysis_policy_reuse policy can greatly improve computation performance of meta data. However, the user need to make sure that the sparsity pattern remains unchanged. If this cannot be assured, rocsparse_analysis_policy_force has to be used.- Note - If the matrix sparsity pattern changes, the gathered information will become invalid. - Note - This function is blocking with respect to the host. - Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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. 
- descr – [in] descriptor of the sparse CSR matrix. 
- csr_val – [in] array of - nnzelements of the sparse CSR matrix.
- csr_row_ptr – [in] array of - m+1elements that point to the start of every row of the sparse CSR matrix.
- csr_col_ind – [in] array of - nnzelements containing the column indices of the sparse CSR matrix.
- info – [out] structure that holds the information collected during the analysis step. 
- analysis – [in] rocsparse_analysis_policy_reuse or rocsparse_analysis_policy_force. 
- solve – [in] rocsparse_solve_policy_auto. 
- temp_buffer – [in] temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mor- nnzis invalid.
- rocsparse_status_invalid_pointer – - descr,- csr_val,- csr_row_ptr,- csr_col_ind,- infoor- temp_bufferpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_not_implemented – - trans!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
 
 
rocsparse_csrilu0()#
- 
rocsparse_status rocsparse_scsrilu0(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, float *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- 
rocsparse_status rocsparse_dcsrilu0(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, double *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- 
rocsparse_status rocsparse_ccsrilu0(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, rocsparse_float_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- 
rocsparse_status rocsparse_zcsrilu0(rocsparse_handle handle, rocsparse_int m, rocsparse_int nnz, const rocsparse_mat_descr descr, rocsparse_double_complex *csr_val, const rocsparse_int *csr_row_ptr, const rocsparse_int *csr_col_ind, rocsparse_mat_info info, rocsparse_solve_policy policy, void *temp_buffer)#
- Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csrilu0computes 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 \]- rocsparse_csrilu0requires a user allocated temporary buffer. Its size is returned by rocsparse_scsrilu0_buffer_size(), rocsparse_dcsrilu0_buffer_size(), rocsparse_ccsrilu0_buffer_size() or rocsparse_zcsrilu0_buffer_size(). Furthermore, analysis meta data is required. It can be obtained by rocsparse_scsrilu0_analysis(), rocsparse_dcsrilu0_analysis(), rocsparse_ccsrilu0_analysis() or rocsparse_zcsrilu0_analysis().- rocsparse_csrilu0reports the first zero pivot (either numerical or structural zero). The zero pivot status can be obtained by calling rocsparse_csrilu0_zero_pivot().- Example
- Consider the sparse \(m \times m\) matrix \(A\), stored in CSR storage format. The following example computes the incomplete LU factorization \(M \approx LU\) and solves the preconditioned system \(My = x\). - // Create rocSPARSE handle rocsparse_handle handle; rocsparse_create_handle(&handle); // Create matrix descriptor for M rocsparse_mat_descr descr_M; rocsparse_create_mat_descr(&descr_M); // Create matrix descriptor for L rocsparse_mat_descr descr_L; rocsparse_create_mat_descr(&descr_L); rocsparse_set_mat_fill_mode(descr_L, rocsparse_fill_mode_lower); rocsparse_set_mat_diag_type(descr_L, rocsparse_diag_type_unit); // Create matrix descriptor for U rocsparse_mat_descr descr_U; rocsparse_create_mat_descr(&descr_U); rocsparse_set_mat_fill_mode(descr_U, rocsparse_fill_mode_upper); rocsparse_set_mat_diag_type(descr_U, rocsparse_diag_type_non_unit); // Create matrix info structure rocsparse_mat_info info; rocsparse_create_mat_info(&info); // Obtain required buffer size size_t buffer_size_M; size_t buffer_size_L; size_t buffer_size_U; rocsparse_dcsrilu0_buffer_size(handle, m, nnz, descr_M, csr_val, csr_row_ptr, csr_col_ind, info, &buffer_size_M); rocsparse_dcsrsv_buffer_size(handle, rocsparse_operation_none, m, nnz, descr_L, csr_val, csr_row_ptr, csr_col_ind, info, &buffer_size_L); rocsparse_dcsrsv_buffer_size(handle, rocsparse_operation_none, m, nnz, descr_U, csr_val, csr_row_ptr, csr_col_ind, info, &buffer_size_U); size_t buffer_size = max(buffer_size_M, max(buffer_size_L, buffer_size_U)); // Allocate temporary buffer void* temp_buffer; hipMalloc(&temp_buffer, buffer_size); // Perform analysis steps, using rocsparse_analysis_policy_reuse to improve // computation performance rocsparse_dcsrilu0_analysis(handle, m, nnz, descr_M, csr_val, csr_row_ptr, csr_col_ind, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); rocsparse_dcsrsv_analysis(handle, rocsparse_operation_none, m, nnz, descr_L, csr_val, csr_row_ptr, csr_col_ind, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); rocsparse_dcsrsv_analysis(handle, rocsparse_operation_none, m, nnz, descr_U, csr_val, csr_row_ptr, csr_col_ind, info, rocsparse_analysis_policy_reuse, rocsparse_solve_policy_auto, temp_buffer); // Check for zero pivot rocsparse_int position; if(rocsparse_status_zero_pivot == rocsparse_csrilu0_zero_pivot(handle, info, &position)) { printf("A has structural zero at A(%d,%d)\n", position, position); } // Compute incomplete LU factorization rocsparse_dcsrilu0(handle, m, nnz, descr_M, csr_val, csr_row_ptr, csr_col_ind, info, rocsparse_solve_policy_auto, temp_buffer); // Check for zero pivot if(rocsparse_status_zero_pivot == rocsparse_csrilu0_zero_pivot(handle, info, &position)) { printf("U has structural and/or numerical zero at U(%d,%d)\n", position, position); } // Solve Lz = x rocsparse_dcsrsv_solve(handle, rocsparse_operation_none, m, nnz, &alpha, descr_L, csr_val, csr_row_ptr, csr_col_ind, info, x, z, rocsparse_solve_policy_auto, temp_buffer); // Solve Uy = z rocsparse_dcsrsv_solve(handle, rocsparse_operation_none, m, nnz, &alpha, descr_U, csr_val, csr_row_ptr, csr_col_ind, info, z, y, rocsparse_solve_policy_auto, temp_buffer); // Clean up hipFree(temp_buffer); rocsparse_destroy_mat_info(info); rocsparse_destroy_mat_descr(descr_M); rocsparse_destroy_mat_descr(descr_L); rocsparse_destroy_mat_descr(descr_U); rocsparse_destroy_handle(handle); 
 - Note - The sparse CSR matrix has to be sorted. This can be achieved by calling rocsparse_csrsort(). - Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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. 
- descr – [in] descriptor of the sparse CSR matrix. 
- csr_val – [inout] array of - nnzelements of the sparse CSR matrix.
- csr_row_ptr – [in] array of - m+1elements that point to the start of every row of the sparse CSR matrix.
- csr_col_ind – [in] array of - nnzelements containing the column indices of the sparse CSR matrix.
- info – [in] structure that holds the information collected during the analysis step. 
- policy – [in] rocsparse_solve_policy_auto. 
- temp_buffer – [in] temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mor- nnzis invalid.
- rocsparse_status_invalid_pointer – - descr,- csr_val,- csr_row_ptror- csr_col_indpointer is invalid.
- rocsparse_status_arch_mismatch – the device is not supported. 
- rocsparse_status_internal_error – an internal error occurred. 
- rocsparse_status_not_implemented – - trans!= rocsparse_operation_none or rocsparse_matrix_type != rocsparse_matrix_type_general.
 
 
rocsparse_csrilu0_clear()#
- 
rocsparse_status rocsparse_csrilu0_clear(rocsparse_handle handle, rocsparse_mat_info info)#
- Incomplete LU factorization with 0 fill-ins and no pivoting using CSR storage format. - rocsparse_csrilu0_cleardeallocates all memory that was allocated by rocsparse_scsrilu0_analysis(), rocsparse_dcsrilu0_analysis(), rocsparse_ccsrilu0_analysis() or rocsparse_zcsrilu0_analysis(). This is especially useful, if memory is an issue and the analysis data is not required for further computation.- Note - Calling - rocsparse_csrilu0_clearis optional. All allocated resources will be cleared, when the opaque rocsparse_mat_info struct is destroyed using rocsparse_destroy_mat_info().- Note - This routine does not support execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- info – [inout] structure that holds the information collected during the analysis step. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_pointer – - infopointer is invalid.
- rocsparse_status_memory_error – the buffer holding the meta data could not be deallocated. 
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_gtsv_buffer_size()#
- 
rocsparse_status rocsparse_sgtsv_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const float *dl, const float *d, const float *du, const float *B, rocsparse_int ldb, size_t *buffer_size)#
- 
rocsparse_status rocsparse_dgtsv_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const double *dl, const double *d, const double *du, const double *B, rocsparse_int ldb, size_t *buffer_size)#
- 
rocsparse_status rocsparse_cgtsv_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const rocsparse_float_complex *dl, const rocsparse_float_complex *d, const rocsparse_float_complex *du, const rocsparse_float_complex *B, rocsparse_int ldb, size_t *buffer_size)#
- 
rocsparse_status rocsparse_zgtsv_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const rocsparse_double_complex *dl, const rocsparse_double_complex *d, const rocsparse_double_complex *du, const rocsparse_double_complex *B, rocsparse_int ldb, size_t *buffer_size)#
- Tridiagonal solver with pivoting. - rocsparse_gtsv_buffer_sizereturns the size of the temporary storage buffer that is required by rocsparse_sgtsv(), rocsparse_dgtsv(), rocsparse_cgtsv() and rocsparse_zgtsv(). The temporary storage buffer must be allocated by the user.- Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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).
- buffer_size – [out] number of bytes of the temporary storage buffer required by rocsparse_sgtsv(), rocsparse_dgtsv(), rocsparse_cgtsv() and rocsparse_zgtsv(). 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - m,- nor- ldbis invalid.
- rocsparse_status_invalid_pointer – - dl,- d,- du,- Bor- buffer_sizepointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_gtsv()#
- 
rocsparse_status rocsparse_sgtsv(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const float *dl, const float *d, const float *du, float *B, rocsparse_int ldb, void *temp_buffer)#
- 
rocsparse_status rocsparse_dgtsv(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const double *dl, const double *d, const double *du, double *B, rocsparse_int ldb, void *temp_buffer)#
- 
rocsparse_status rocsparse_cgtsv(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const rocsparse_float_complex *dl, const rocsparse_float_complex *d, const rocsparse_float_complex *du, rocsparse_float_complex *B, rocsparse_int ldb, void *temp_buffer)#
- 
rocsparse_status rocsparse_zgtsv(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const rocsparse_double_complex *dl, const rocsparse_double_complex *d, const rocsparse_double_complex *du, rocsparse_double_complex *B, rocsparse_int ldb, void *temp_buffer)#
- Tridiagonal solver with pivoting. - rocsparse_gtsvsolves 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. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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).
- temp_buffer – [in] temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - m,- nor- ldbis invalid.
- rocsparse_status_invalid_pointer – - dl,- d,- du,- Bor- temp_bufferpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_gtsv_no_pivot_buffer_size()#
- 
rocsparse_status rocsparse_sgtsv_no_pivot_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const float *dl, const float *d, const float *du, const float *B, rocsparse_int ldb, size_t *buffer_size)#
- 
rocsparse_status rocsparse_dgtsv_no_pivot_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const double *dl, const double *d, const double *du, const double *B, rocsparse_int ldb, size_t *buffer_size)#
- 
rocsparse_status rocsparse_cgtsv_no_pivot_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const rocsparse_float_complex *dl, const rocsparse_float_complex *d, const rocsparse_float_complex *du, const rocsparse_float_complex *B, rocsparse_int ldb, size_t *buffer_size)#
- 
rocsparse_status rocsparse_zgtsv_no_pivot_buffer_size(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const rocsparse_double_complex *dl, const rocsparse_double_complex *d, const rocsparse_double_complex *du, const rocsparse_double_complex *B, rocsparse_int ldb, size_t *buffer_size)#
- Tridiagonal solver (no pivoting) - rocsparse_gtsv_no_pivot_buffer_sizereturns the size of the temporary storage buffer that is required by rocsparse_sgtsv_no_pivot(), rocsparse_dgtsv_no_pivot(), rocsparse_cgtsv_no_pivot() and rocsparse_zgtsv_no_pivot(). The temporary storage buffer must be allocated by the user.- Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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).
- buffer_size – [out] number of bytes of the temporary storage buffer required by rocsparse_sgtsv_no_pivot(), rocsparse_dgtsv_no_pivot(), rocsparse_cgtsv_no_pivot() and rocsparse_zgtsv_no_pivot(). 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - m,- nor- ldbis invalid.
- rocsparse_status_invalid_pointer – - dl,- d,- du,- Bor- buffer_sizepointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_gtsv_no_pivot()#
- 
rocsparse_status rocsparse_sgtsv_no_pivot(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const float *dl, const float *d, const float *du, float *B, rocsparse_int ldb, void *temp_buffer)#
- 
rocsparse_status rocsparse_dgtsv_no_pivot(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const double *dl, const double *d, const double *du, double *B, rocsparse_int ldb, void *temp_buffer)#
- 
rocsparse_status rocsparse_cgtsv_no_pivot(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const rocsparse_float_complex *dl, const rocsparse_float_complex *d, const rocsparse_float_complex *du, rocsparse_float_complex *B, rocsparse_int ldb, void *temp_buffer)#
- 
rocsparse_status rocsparse_zgtsv_no_pivot(rocsparse_handle handle, rocsparse_int m, rocsparse_int n, const rocsparse_double_complex *dl, const rocsparse_double_complex *d, const rocsparse_double_complex *du, rocsparse_double_complex *B, rocsparse_int ldb, void *temp_buffer)#
- Tridiagonal solver (no pivoting) - rocsparse_gtsv_no_pivotsolves 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. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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).
- temp_buffer – [in] temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - m,- nor- ldbis invalid.
- rocsparse_status_invalid_pointer – - dl,- d,- du,- Bor- temp_bufferpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_gtsv_no_pivot_strided_batch_buffer_size()#
- 
rocsparse_status rocsparse_sgtsv_no_pivot_strided_batch_buffer_size(rocsparse_handle handle, rocsparse_int m, const float *dl, const float *d, const float *du, const float *x, rocsparse_int batch_count, rocsparse_int batch_stride, size_t *buffer_size)#
- 
rocsparse_status rocsparse_dgtsv_no_pivot_strided_batch_buffer_size(rocsparse_handle handle, rocsparse_int m, const double *dl, const double *d, const double *du, const double *x, rocsparse_int batch_count, rocsparse_int batch_stride, size_t *buffer_size)#
- 
rocsparse_status rocsparse_cgtsv_no_pivot_strided_batch_buffer_size(rocsparse_handle handle, rocsparse_int m, const rocsparse_float_complex *dl, const rocsparse_float_complex *d, const rocsparse_float_complex *du, const rocsparse_float_complex *x, rocsparse_int batch_count, rocsparse_int batch_stride, size_t *buffer_size)#
- 
rocsparse_status rocsparse_zgtsv_no_pivot_strided_batch_buffer_size(rocsparse_handle handle, rocsparse_int m, const rocsparse_double_complex *dl, const rocsparse_double_complex *d, const rocsparse_double_complex *du, const rocsparse_double_complex *x, rocsparse_int batch_count, rocsparse_int batch_stride, size_t *buffer_size)#
- Strided Batch tridiagonal solver (no pivoting) - rocsparse_gtsv_no_pivot_strided_batch_buffer_sizereturns the size of the temporary storage buffer that is required by rocsparse_sgtsv_no_pivot_strided_batch(), rocsparse_dgtsv_no_pivot_strided_batch(), rocsparse_cgtsv_no_pivot_strided_batch() and rocsparse_zgtsv_no_pivot_strided_batch(). The temporary storage buffer must be allocated by the user.- Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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+batch_stride*i.
- d – [in] main diagonal of tri-diagonal system where the ith system diagonal starts at - d+batch_stride*i.
- du – [in] upper diagonal of tri-diagonal system where the ith system upper diagonal starts at - du+batch_stride*i.
- x – [inout] Dense array of righthand-sides where the ith righthand-side starts at - x+batch_stride*i.
- batch_count – [in] The number of systems to solve. 
- batch_stride – [in] The number of elements that separate each system. Must satisfy - batch_stride>= m.
- buffer_size – [out] number of bytes of the temporary storage buffer required by rocsparse_sgtsv_no_pivot_strided_batch(), rocsparse_dgtsv_no_pivot_strided_batch(), rocsparse_cgtsv_no_pivot_strided_batch() and rocsparse_zgtsv_no_pivot_strided_batch(). 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - m,- batch_countor- batch_strideis invalid.
- rocsparse_status_invalid_pointer – - dl,- d,- du,- xor- buffer_sizepointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_gtsv_no_pivot_strided_batch()#
- 
rocsparse_status rocsparse_sgtsv_no_pivot_strided_batch(rocsparse_handle handle, rocsparse_int m, const float *dl, const float *d, const float *du, float *x, rocsparse_int batch_count, rocsparse_int batch_stride, void *temp_buffer)#
- 
rocsparse_status rocsparse_dgtsv_no_pivot_strided_batch(rocsparse_handle handle, rocsparse_int m, const double *dl, const double *d, const double *du, double *x, rocsparse_int batch_count, rocsparse_int batch_stride, void *temp_buffer)#
- 
rocsparse_status rocsparse_cgtsv_no_pivot_strided_batch(rocsparse_handle handle, rocsparse_int m, const rocsparse_float_complex *dl, const rocsparse_float_complex *d, const rocsparse_float_complex *du, rocsparse_float_complex *x, rocsparse_int batch_count, rocsparse_int batch_stride, void *temp_buffer)#
- 
rocsparse_status rocsparse_zgtsv_no_pivot_strided_batch(rocsparse_handle handle, rocsparse_int m, const rocsparse_double_complex *dl, const rocsparse_double_complex *d, const rocsparse_double_complex *du, rocsparse_double_complex *x, rocsparse_int batch_count, rocsparse_int batch_stride, void *temp_buffer)#
- Strided Batch tridiagonal solver (no pivoting) - rocsparse_gtsv_no_pivot_strided_batchsolves 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. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse 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+batch_stride*i.
- batch_count – [in] The number of systems to solve. 
- batch_stride – [in] The number of elements that separate each system. Must satisfy - batch_stride>= m.
- temp_buffer – [in] temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - m,- batch_countor- batch_strideis invalid.
- rocsparse_status_invalid_pointer – - dl,- d,- du,- xor- temp_bufferpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_gtsv_interleaved_batch_buffer_size()#
- 
rocsparse_status rocsparse_sgtsv_interleaved_batch_buffer_size(rocsparse_handle handle, rocsparse_gtsv_interleaved_alg alg, rocsparse_int m, const float *dl, const float *d, const float *du, const float *x, rocsparse_int batch_count, rocsparse_int batch_stride, size_t *buffer_size)#
- 
rocsparse_status rocsparse_dgtsv_interleaved_batch_buffer_size(rocsparse_handle handle, rocsparse_gtsv_interleaved_alg alg, rocsparse_int m, const double *dl, const double *d, const double *du, const double *x, rocsparse_int batch_count, rocsparse_int batch_stride, size_t *buffer_size)#
- 
rocsparse_status rocsparse_cgtsv_interleaved_batch_buffer_size(rocsparse_handle handle, rocsparse_gtsv_interleaved_alg alg, rocsparse_int m, const rocsparse_float_complex *dl, const rocsparse_float_complex *d, const rocsparse_float_complex *du, const rocsparse_float_complex *x, rocsparse_int batch_count, rocsparse_int batch_stride, size_t *buffer_size)#
- 
rocsparse_status rocsparse_zgtsv_interleaved_batch_buffer_size(rocsparse_handle handle, rocsparse_gtsv_interleaved_alg alg, rocsparse_int m, const rocsparse_double_complex *dl, const rocsparse_double_complex *d, const rocsparse_double_complex *du, const rocsparse_double_complex *x, rocsparse_int batch_count, rocsparse_int batch_stride, size_t *buffer_size)#
- Interleaved Batch tridiagonal solver. - rocsparse_gtsv_interleaved_batch_buffer_sizereturns the size of the temporary storage buffer that is required by rocsparse_sgtsv_interleaved_batch(), rocsparse_dgtsv_interleaved_batch(), rocsparse_cgtsv_interleaved_batch() and rocsparse_zgtsv_interleaved_batch(). The temporary storage buffer must be allocated by the user.- Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- alg – [in] Algorithm to use when solving tridiagonal systems. Options are thomas ( - rocsparse_gtsv_interleaved_thomas), LU (- rocsparse_gtsv_interleaved_lu), or QR (- rocsparse_gtsv_interleaved_qr). Passing- rocsparse_gtsv_interleaved_defaultdefaults the algorithm to use QR. 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 - batch_strideby- m.
- batch_count – [in] The number of systems to solve. 
- batch_stride – [in] The number of elements that separate consecutive elements in a system. Must satisfy - batch_stride>= batch_count.
- buffer_size – [out] number of bytes of the temporary storage buffer required by rocsparse_sgtsv_interleaved_batch(), rocsparse_dgtsv_interleaved_batch(), rocsparse_cgtsv_interleaved_batch() and rocsparse_zgtsv_interleaved_batch(). 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - m,- batch_count,- batch_strideis invalid.
- rocsparse_status_invalid_pointer – - dl,- d,- du,- xor- buffer_sizepointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_gtsv_interleaved_batch()#
- 
rocsparse_status rocsparse_sgtsv_interleaved_batch(rocsparse_handle handle, rocsparse_gtsv_interleaved_alg alg, rocsparse_int m, float *dl, float *d, float *du, float *x, rocsparse_int batch_count, rocsparse_int batch_stride, void *temp_buffer)#
- 
rocsparse_status rocsparse_dgtsv_interleaved_batch(rocsparse_handle handle, rocsparse_gtsv_interleaved_alg alg, rocsparse_int m, double *dl, double *d, double *du, double *x, rocsparse_int batch_count, rocsparse_int batch_stride, void *temp_buffer)#
- 
rocsparse_status rocsparse_cgtsv_interleaved_batch(rocsparse_handle handle, rocsparse_gtsv_interleaved_alg alg, rocsparse_int m, rocsparse_float_complex *dl, rocsparse_float_complex *d, rocsparse_float_complex *du, rocsparse_float_complex *x, rocsparse_int batch_count, rocsparse_int batch_stride, void *temp_buffer)#
- 
rocsparse_status rocsparse_zgtsv_interleaved_batch(rocsparse_handle handle, rocsparse_gtsv_interleaved_alg alg, rocsparse_int m, rocsparse_double_complex *dl, rocsparse_double_complex *d, rocsparse_double_complex *du, rocsparse_double_complex *x, rocsparse_int batch_count, rocsparse_int batch_stride, void *temp_buffer)#
- Interleaved Batch tridiagonal solver. - rocsparse_gtsv_interleaved_batchsolves a batched tridiagonal linear system. The routine requires a temporary storage buffer that must be allocated by the user. The size of this buffer can be determined by first calling- rocsparse_gtsv_interleaved_batch_buffer_size. The user can specify different algorithms for- rocsparse_gtsv_interleaved_batchto use. Options are thomas (- rocsparse_gtsv_interleaved_thomas), LU (- rocsparse_gtsv_interleaved_lu), or QR (- rocsparse_gtsv_interleaved_qr). Passing- rocsparse_gtsv_interleaved_defaultdefaults the algorithm to use QR.- Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- alg – [in] Algorithm to use when solving tridiagonal systems. Options are thomas ( - rocsparse_gtsv_interleaved_thomas), LU (- rocsparse_gtsv_interleaved_lu), or QR (- rocsparse_gtsv_interleaved_qr). Passing- rocsparse_gtsv_interleaved_defaultdefaults the algorithm to use QR. 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 - batch_strideby- m.
- batch_count – [in] The number of systems to solve. 
- batch_stride – [in] The number of elements that separate consecutive elements in a system. Must satisfy - batch_stride>= batch_count.
- temp_buffer – [in] temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - mor- batch_countor- batch_strideis invalid.
- rocsparse_status_invalid_pointer – - dl,- d,- du,- xor- temp_bufferpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_gpsv_interleaved_batch_buffer_size()#
- 
rocsparse_status rocsparse_sgpsv_interleaved_batch_buffer_size(rocsparse_handle handle, rocsparse_gpsv_interleaved_alg alg, rocsparse_int m, const float *ds, const float *dl, const float *d, const float *du, const float *dw, const float *x, rocsparse_int batch_count, rocsparse_int batch_stride, size_t *buffer_size)#
- 
rocsparse_status rocsparse_dgpsv_interleaved_batch_buffer_size(rocsparse_handle handle, rocsparse_gpsv_interleaved_alg alg, rocsparse_int m, const double *ds, const double *dl, const double *d, const double *du, const double *dw, const double *x, rocsparse_int batch_count, rocsparse_int batch_stride, size_t *buffer_size)#
- 
rocsparse_status rocsparse_cgpsv_interleaved_batch_buffer_size(rocsparse_handle handle, rocsparse_gpsv_interleaved_alg alg, rocsparse_int m, const rocsparse_float_complex *ds, const rocsparse_float_complex *dl, const rocsparse_float_complex *d, const rocsparse_float_complex *du, const rocsparse_float_complex *dw, const rocsparse_float_complex *x, rocsparse_int batch_count, rocsparse_int batch_stride, size_t *buffer_size)#
- 
rocsparse_status rocsparse_zgpsv_interleaved_batch_buffer_size(rocsparse_handle handle, rocsparse_gpsv_interleaved_alg alg, rocsparse_int m, const rocsparse_double_complex *ds, const rocsparse_double_complex *dl, const rocsparse_double_complex *d, const rocsparse_double_complex *du, const rocsparse_double_complex *dw, const rocsparse_double_complex *x, rocsparse_int batch_count, rocsparse_int batch_stride, size_t *buffer_size)#
- Batched Pentadiagonal solver. - rocsparse_gpsv_interleaved_batch_buffer_sizecalculates the required buffer size for rocsparse_gpsv_interleaved_batch(). It is the users responsibility to allocate this buffer.- Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- alg – [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 - batch_strideby- m.
- batch_count – [in] The number of systems to solve. 
- batch_stride – [in] The number of elements that separate consecutive elements in a system. Must satisfy - batch_stride>= batch_count.
- buffer_size – [out] Number of bytes of the temporary storage buffer required. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - m,- alg,- batch_countor- batch_strideis invalid.
- rocsparse_status_invalid_pointer – - ds,- dl,- d,- du,- dw,- xor- temp_bufferpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred. 
 
 
rocsparse_gpsv_interleaved_batch()#
- 
rocsparse_status rocsparse_sgpsv_interleaved_batch(rocsparse_handle handle, rocsparse_gpsv_interleaved_alg alg, rocsparse_int m, float *ds, float *dl, float *d, float *du, float *dw, float *x, rocsparse_int batch_count, rocsparse_int batch_stride, void *temp_buffer)#
- 
rocsparse_status rocsparse_dgpsv_interleaved_batch(rocsparse_handle handle, rocsparse_gpsv_interleaved_alg alg, rocsparse_int m, double *ds, double *dl, double *d, double *du, double *dw, double *x, rocsparse_int batch_count, rocsparse_int batch_stride, void *temp_buffer)#
- 
rocsparse_status rocsparse_cgpsv_interleaved_batch(rocsparse_handle handle, rocsparse_gpsv_interleaved_alg alg, rocsparse_int m, rocsparse_float_complex *ds, rocsparse_float_complex *dl, rocsparse_float_complex *d, rocsparse_float_complex *du, rocsparse_float_complex *dw, rocsparse_float_complex *x, rocsparse_int batch_count, rocsparse_int batch_stride, void *temp_buffer)#
- 
rocsparse_status rocsparse_zgpsv_interleaved_batch(rocsparse_handle handle, rocsparse_gpsv_interleaved_alg alg, rocsparse_int m, rocsparse_double_complex *ds, rocsparse_double_complex *dl, rocsparse_double_complex *d, rocsparse_double_complex *du, rocsparse_double_complex *dw, rocsparse_double_complex *x, rocsparse_int batch_count, rocsparse_int batch_stride, void *temp_buffer)#
- Batched Pentadiagonal solver. - rocsparse_gpsv_interleaved_batchsolves a batch of pentadiagonal linear systems. The coefficient matrix of each pentadiagonal linear system is defined by five vectors for the lower part (ds, dl), main diagonal (d) and upper part (du, dw).- The function requires a temporary buffer. The size of the required buffer is returned by rocsparse_gpsv_interleaved_batch_buffer_size(). - Note - This function is non blocking and executed asynchronously with respect to the host. It may return before the actual computation has finished. - Note - The routine is numerically stable because it uses QR to solve the linear systems. - Note - m need to be at least 3, to be a valid pentadiagonal matrix. - Note - This routine supports execution in a hipGraph context. - Parameters:
- handle – [in] handle to the rocsparse library context queue. 
- alg – [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 - batch_strideby- m.
- batch_count – [in] The number of systems to solve. 
- batch_stride – [in] The number of elements that separate consecutive elements in a system. Must satisfy - batch_stride>= batch_count.
- temp_buffer – [in] Temporary storage buffer allocated by the user. 
 
- Return values:
- rocsparse_status_success – the operation completed successfully. 
- rocsparse_status_invalid_handle – the library context was not initialized. 
- rocsparse_status_invalid_size – - m,- alg,- batch_countor- batch_strideis invalid.
- rocsparse_status_invalid_pointer – - ds,- dl,- d,- du,- dw,- xor- temp_bufferpointer is invalid.
- rocsparse_status_internal_error – an internal error occurred.