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

?syev

Computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix.

Syntax

lapack_int LAPACKE_ssyev (int matrix_layout, char jobz, char uplo, lapack_int n, float* a, lapack_int lda, float* w);

lapack_int LAPACKE_dsyev (int matrix_layout, char jobz, char uplo, lapack_int n, double* a, lapack_int lda, double* w);

Include Files

  • mkl.h

Description

The routine computes all eigenvalues and, optionally, eigenvectors of a real symmetric matrix A.

Note that for most cases of real symmetric eigenvalue problems the default choice should be syevr function as its underlying algorithm is faster and uses less workspace.

Input Parameters

matrix_layout

Specifies whether matrix storage layout is row major (LAPACK_ROW_MAJOR) or column major (LAPACK_COL_MAJOR).

jobz

Must be 'N' or 'V'.

If jobz = 'N', then only eigenvalues are computed.

If jobz = 'V', then eigenvalues and eigenvectors are computed.

uplo

Must be 'U' or 'L'.

If uplo = 'U', a stores the upper triangular part of A.

If uplo = 'L', a stores the lower triangular part of A.

n

The order of the matrix A (n 0).

a

a (size max(1, lda*n)) is an array containing either upper or lower triangular part of the symmetric matrix A, as specified by uplo.

lda

The leading dimension of the array a.

Must be at least max(1, n).

Output Parameters

a

On exit, if jobz = 'V', then if info = 0, array a contains the orthonormal eigenvectors of the matrix A.

If jobz = 'N', then on exit the lower triangle

(if uplo = 'L') or the upper triangle (if uplo = 'U') of A, including the diagonal, is overwritten.

w

Array, size at least max(1, n).

If info = 0, contains the eigenvalues of the matrix A in ascending order.

Return Values

This function returns a value info.

If info=0, the execution is successful.

If info = -i, the i-th parameter had an illegal value.

If info = i, then the algorithm failed to converge; i indicates the number of elements of an intermediate tridiagonal form which did not converge to zero.