Developer Reference for Intel® oneAPI Math Kernel Library for C

ID 766684
Date 11/07/2023
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_sypr

Computes the symmetric product of three sparse matrices and stores the result in a newly allocated sparse matrix.

Syntax

sparse_status_t mkl_sparse_sypr (const sparse_operation_t operation , const sparse_matrix_t A, const sparse_matrix_t B, const struct matrix_descr B, sparse_matrix_t *C, const sparse_request_t request);

Include Files

  • mkl_spblas.h

Description

The mkl_sparse_sypr routine performs a multiplication of three sparse matrices that results in a symmetric or Hermitian matrix, C.

C:=A*B*opA(A)
or
C:=opA(A)*B*A
depending on the matrix modifier operation.

Here, A, B, and C are sparse matrices, where A has a general structure while B and C are symmetric (for real data types) or Hermitian (for complex data types) matrices. opA is the transpose (real data types) or conjugate transpose (complex data types) operator.

NOTE:

This routine is not supported for sparse matrices in COO or CSC formats. This routine supports only CSR and BSR formats. In addition, it supports only the sorted CSR and sorted BSR formats for the input matrix. If the data is unsorted, call the mkl_sparse_order routine before either mkl_sparse_sypr or mkl_sparse_?_syprd.

Input Parameters

operation

Specifies operation on the input sparse matrices.

SPARSE_OPERATION_NON_TRANSPOSE

Non-transpose case.

C:=A*B*(AT) for real precision

C:=A*B*(AH) for complex precision.

SPARSE_OPERATION_TRANSPOSE

Transpose case. This is not supported for complex matrices.

C:=(AT)*B*A

SPARSE_OPERATION_CONJUGATE_TRANSPOSE

Conjugate transpose case. This is not supported for real matrices.

C:=(AH)*B*A

A

Handle which contains the sparse matrix A.

B

Handle which contains the sparse matrix B.

descrB

Structure specifying properties of the sparse matrix.

sparse_matrix_type_t type specifies the type of a sparse matrix

SPARSE_MATRIX_TYPE_SYMMETRIC

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

SPARSE_MATRIX_TYPE_HERMITIAN

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

sparse_fill_mode_t mode specifies the triangular matrix part.

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_t diag specifies the type of diagonal.

SPARSE_DIAG_NON_UNIT

Diagonal elements cannot be equal to one.

NOTE:

This routine also supports C=AAT,H with these parameters:

descrB.type=SPARSE_MATRIX_TYPE_DIAGONAL

descrB.diag=SPARSE_DIAG_UNIT

In this case, you do not need to allocate structure B. Use the routine as a 2-stage version of mkl_sparse_syrk.

request

Use this routine to specify if the computations should be performed in a single step or using the two-stage algorithm. See Two-stage Algorithm for Inspector-executor Sparse BLAS Routines for more information.

SPARSE_STAGE_NNZ_COUNT Only rowIndex (BSR/CSR format) or colIndex (CSC format) array of the matrix is computed internally. The computation can be extracted to measure the memory required for full operation.
SPARSE_STAGE_FINALIZE_MULT_NO_VAL Finalize computations of the matrix structure (values will not be computed). Use only after the call with SPARSE_STAGE_NNZ_COUNT parameter.
SPARSE_STAGE_FINALIZE_MULT Finalize computation. Can be used after the call with the SPARSE_STAGE_NNZ_COUNT or SPARSE_STAGE_FINALIZE_MULT_NO_VAL. Can also be used when the matrix structure remains unchanged and only values of the resulting matrix C need to be recomputed.
SPARSE_STAGE_FULL_MULT_NO_VAL Perform computations of the matrix structure.
SPARSE_STAGE_FULL_MULT Perform the entire computation in a single step.

Output Parameters

C

Handle which contains the resulting sparse matrix. Only the upper-triangular part of the matrix is computed.

Return Values

The function returns a value indicating whether the operation was successful, or the reason why it failed.

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.