Visible to Intel only — GUID: GUID-6E6360FD-4FE6-41C9-9092-BC7A67ED7967
Visible to Intel only — GUID: GUID-6E6360FD-4FE6-41C9-9092-BC7A67ED7967
?hbtrd
Reduces a complex Hermitian band matrix to tridiagonal form.
lapack_int LAPACKE_chbtrd( int matrix_layout, char vect, char uplo, lapack_int n, lapack_int kd, lapack_complex_float* ab, lapack_int ldab, float* d, float* e, lapack_complex_float* q, lapack_int ldq );
lapack_int LAPACKE_zhbtrd( int matrix_layout, char vect, char uplo, lapack_int n, lapack_int kd, lapack_complex_double* ab, lapack_int ldab, double* d, double* e, lapack_complex_double* q, lapack_int ldq );
- mkl.h
The routine reduces a complex Hermitian band matrix A to symmetric tridiagonal form T by a unitary similarity transformation: A = Q*T*QH. The unitary matrix Q is determined as a product of Givens rotations.
If required, the routine can also form the matrix Q explicitly.
- matrix_layout
-
Specifies whether matrix storage layout is row major (LAPACK_ROW_MAJOR) or column major (LAPACK_COL_MAJOR).
- vect
-
Must be 'V', 'N', or 'U'.
If vect = 'V', the routine returns the explicit matrix Q.
If vect = 'N', the routine does not return Q.
If vect = 'U', the routine updates matrix X by forming Q*X.
- uplo
-
Must be 'U' or 'L'.
If uplo = 'U', ab stores the upper triangular part of A.
If uplo = 'L', ab stores the lower triangular part of A.
- n
-
The order of the matrix A (n≥ 0).
- kd
-
The number of super- or sub-diagonals in A
(kd≥ 0).
- ab
-
ab(size at least max(1, ldab*n) for column major layout and at least max(1, ldab*(kd+ 1)) for row major layout) is an array containing either upper or lower triangular part of the matrix A (as specified by uplo) in band storage format.
- q
-
q (size max(1, ldq*n)) is an array.
If vect = 'U', the q array must contain an n-by-n matrix X.
If vect = 'N' or 'V', the q parameter need not be set.'
- ldab
-
The leading dimension of ab; at least kd+1 for column major layout and n for row major layout.
- ldq
-
The leading dimension of q. Constraints:
ldq≥ max(1, n) if vect = 'V' or 'U';
ldq≥ 1 if vect = 'N'.
- ab
-
On exit, the diagonal elements of the array ab are overwritten by the diagonal elements of the tridiagonal matrix T. If kd > 0, the elements on the first superdiagonal (if uplo = 'U') or the first subdiagonal (if uplo = 'L') are ovewritten by the off-diagonal elements of T. The rest of ab is overwritten by values generated during the reduction.
- d, e
-
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).
- q
-
If vect = 'N', q is not referenced.
If vect = 'V', q contains the n-by-n matrix Q.
If vect = 'U', q contains the product X* Q.
This function returns a value info.
If info=0, the execution is successful.
If info = -i, the i-th parameter had an illegal value.
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 computed matrix Q differs from an exactly unitary matrix by a matrix E such that ||E||2 = O(ε).
The total number of floating-point operations is approximately 20n2*kd if vect = 'N', with 10n3*(kd-1)/kd additional operations if vect = 'V'.
The real counterpart of this routine is sbtrd.