Developer Reference for Intel® oneAPI Math Kernel Library for C

ID 766684
Date 12/16/2022
Public

A newer version of this document is available. Customers should click here to go to the newest version.

Document Table of Contents

mkl_sparse_?_lu_smoother

Computes an action of a preconditioner which corresponds to the approximate matrix decomposition for the system (see description below).

Syntax

sparse_status_t mkl_sparse_s_lu_smoother (const sparse_operation_t op, const sparse_matrix_t A, const struct matrix descr descr, const float *diag, const float *approx_diag_inverse, float *x, const float *b);

sparse_status_t mkl_sparse_d_lu_smoother (const sparse_operation_t op, const sparse_matrix_t A, const struct matrix descr descr, const double *diag, const double *approx_diag_inverse, double *x, const double *b);

sparse_status_t mkl_sparse_c_lu_smoother (const sparse_operation_t op, const sparse_matrix_t A, const struct matrix descr descr, const MKL_COMPLEX8 *diag, const MKL_COMPLEX8 *approx_diag_inverse, MKL_COMPLEX8 *x, const MKL_COMPLEX8 *b);

sparse_status_t mkl_sparse_z_lu_smoother (const sparse_operation_t op, const sparse_matrix_t A, const struct matrix descr descr, const MKL_COMPLEX16 *diag, const MKL_COMPLEX16 *approx_diag_inverse, MKL_COMPLEX16 *x, const MKL_COMPLEX16 *b);

Include Files
  • mkl_spblas.h
Description

This routine computes an update for an iterative solution x of the system Ax=b by means of applying one iteration of an approximate preconditioner which is based on the following approximation:

, where E is an approximate inverse of the diagonal (using exact inverse will result in Gauss-Seidel preconditioner), L and U are lower/upper triangular parts of A, D is the diagonal (block diagonal in case of BSR format) of A.

The mkl_sparse_?_lu_smoother routine performs these operations:

r = b - A*x    /* 1. Computes the residual */
(L + D)*E*(U + D)*dx = r    /* 2. Finds the update dx by solving the system */
y = x + dx    /* 3. Performs an update */

This is also equal to the Symmetric Gauss-Seidel operation in the case of a CSR format and 1x1 diagonal blocks:

(L + D)*x^1 = b - U*x  /* Lower solve for intermediate x^1 */
(U + D)*x = b - L*x^1  /* Upper solve */
NOTE:

This routine is supported only for non-transpose operation, real data types, and CSR/BSR sparse formats. In a BSR format, both diagonal values and approximate diagonal inverse arrays should be passed explicitly. For CSR format, diagonal values should be passed explicitly.

Input Parameters
operation

Specifies the operation performed on matrix A.

SPARSE_OPERATION_NON_TRANSPOSE, op(A) := A

NOTE:

Transpose and conjugate transpose (SPARSE_OPERATION_TRANSPOSE and SPARSE_OPERATION_CONJUGATE_TRANSPOSE) are not supported.

Non-transpose, op(A)= A.

A

Handle which contains the sparse matrix A.

descr

Structure specifying sparse matrix properties.

sparse_matrix_type_ttype - Specifies the type of a sparse matrix:

SPARSE_MATRIX_TYPE_GENERAL

The matrix is processed as is.

SPARSE_MATRIX_TYPE_SYMMETRIC

The matrix is symmetric (only the requested triangle is processed).

SPARSE_MATRIX_TYPE_HERMITIAN

The matrix is Hermitian (only the requested triangle is processed).

SPARSE_MATRIX_TYPE_TRIANGULAR

The matrix is triangular (only the requested triangle is processed).

SPARSE_MATRIX_TYPE_DIAGONAL

The matrix is diagonal (only diagonal elements are processed).

SPARSE_MATRIX_TYPE_BLOCK_TRIANGULAR

The matrix is block-triangular (only the requested triangle is processed). Applies to BSR format only.

SPARSE_MATRIX_TYPE_BLOCK_DIAGONAL

The matrix is block-diagonal (only diagonal blocks are processed). Applies to BSR format only.

sparse_fill_mode_tmode - Specifies the triangular matrix part for symmetric, Hermitian, triangular, and block-triangular matrices:

SPARSE_FILL_MODE_LOWER

The lower triangular matrix part is processed.

SPARSE_FILL_MODE_UPPER

The upper triangular matrix part is processed.

sparse_diag_type_tdiag - Specifies the diagonal type for non-general matrices:

SPARSE_DIAG_NON_UNIT

Diagonal elements might not be equal to one.

SPARSE_DIAG_UNIT

Diagonal elements are equal to one.

NOTE:

Only SPARSE_MATRIX_TYPE_GENERAL is supported.

diag

Array of size at least m, where m is the number of rows (or nrows * block_size * block_size in case of BSR format) of matrix A.

The array diag must contain the diagonal values of matrix A.

approx_diag_inverse

Array of size at least m, where m is the number of rows (or the number of rows * block_size * block_size in case of BSR format) of matrix A.

The array approx_diag_inverse will be used as E, approximate inverse of the diagonal of the matrix A.

x

Array of size at least k, where k is the number of columns (or columns * block_size in case of BSR format) of matrix A.

On entry, the array x must contain the input vector.

b

Array of size at least m, where m is the number of rows ( or rows * block_size in case of BSR format ) of matrix A. The array b must contain the values of the right-hand side of the system.

Output Parameters
x

Overwritten by the computed vector y.

Return Values

The function returns a value indicating whether the operation was successful or not, and why.

SPARSE_STATUS_SUCCESS

The operation was successful.

SPARSE_STATUS_NOT_INITIALIZED

The routine encountered an empty handle or matrix array.

SPARSE_STATUS_ALLOC_FAILED

Internal memory allocation failed.

SPARSE_STATUS_INVALID_VALUE

The input parameters contain an invalid value.

SPARSE_STATUS_EXECUTION_FAILED

Execution failed.

SPARSE_STATUS_INTERNAL_ERROR

An error in algorithm implementation occurred.

SPARSE_STATUS_NOT_SUPPORTED

The requested operation is not supported.