Visible to Intel only — GUID: GUID-3C4516C4-705C-4362-A3B9-C1C5FA0966EB
Visible to Intel only — GUID: GUID-3C4516C4-705C-4362-A3B9-C1C5FA0966EB
?sptrd
Reduces a real symmetric matrix to tridiagonal form using packed storage.
Syntax
lapack_int LAPACKE_ssptrd (int matrix_layout, char uplo, lapack_int n, float* ap, float* d, float* e, float* tau);
lapack_int LAPACKE_dsptrd (int matrix_layout, char uplo, lapack_int n, double* ap, double* d, double* e, double* tau);
Include Files
- mkl.h
Description
The routine reduces a packed real symmetric matrix A to symmetric tridiagonal form T by an orthogonal similarity transformation: A = Q*T*QT. The orthogonal matrix Q is not formed explicitly but is represented as a product of n-1 elementary reflectors. Routines are provided for working with Q in this representation. See Application Notes below for details.
Input Parameters
- matrix_layout
-
Specifies whether matrix storage layout is row major (LAPACK_ROW_MAJOR) or column major (LAPACK_COL_MAJOR).
- uplo
-
Must be 'U' or 'L'.
If uplo = 'U', ap stores the packed upper triangle of A.
If uplo = 'L', ap stores the packed lower triangle of A.
- n
-
The order of the matrix A (n≥ 0).
- ap
-
Array, size at least max(1, n(n+1)/2). Contains either upper or lower triangle of A (as specified by uplo) in the packed form described in Matrix Storage Schemes.
Output Parameters
- ap
-
Overwritten by the tridiagonal matrix T and details of the orthogonal matrix Q, as specified by uplo.
- d, e, tau
-
Arrays:
d contains the diagonal elements of the matrix T.
The dimension of d must be at least max(1, n).
e contains the off-diagonal elements of T.
The dimension of e must be at least max(1, n-1).
tau Stores (n-1) scalars that define elementary reflectors in decomposition of the matrix Q in a product of n-1 reflectors.
The dimension of tau must be at least max(1, n-1).
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.
Application Notes
The matrix Q is represented as a product of n-1 elementary reflectors, as follows :
If uplo = 'U', Q = H(n-1) ... H(2)H(1)
Each H(i) has the form
H(i) = I - tau*v*vT
where tau is a real scalar and v is a real vector with v(i+1:n) = 0 and v(i) = 1.
On exit, tau is stored in tau[i - 1], and v(1:i-1) is stored in AP, overwriting A(1:i-1, i+1).
If uplo = 'L', Q = H(1)H(2) ... H(n-1)
Each H(i) has the form
H(i) = I - tau*v*vT
where tau is a real scalar and v is a real vector with v(1:i) = 0 and v(i+1) = 1.
On exit, tau is stored in tau[i - 1], and v(i+2:n) is stored in AP, overwriting A(i+2:n, i).
The computed matrix T is exactly similar to a matrix A+E, where ||E||2 = c(n)*ε*||A||2, c(n) is a modestly increasing function of n, and ε is the machine precision. The approximate number of floating-point operations is (4/3)n3.
After calling this routine, you can call the following:
The complex counterpart of this routine is hptrd.