hetrd#
Reduces a complex Hermitian matrix to tridiagonal form.
Description
hetrd supports the following precisions.
Routine name
T
chetrd
std::complex<float>
zhetrd
std::complex<double>
The routine reduces a complex Hermitian matrix \(A\) to symmetric tridiagonal form \(T\) by a unitary similarity transformation: \(A = QTQ^H\). The unitary matrix \(Q\) is not formed explicitly but is represented as a product of \(n-1\) elementary reflectors. Routines are provided to work with \(Q\) in this representation.
hetrd (Buffer Version)#
Syntax
namespace oneapi::mkl::lapack {
void hetrd(sycl::queue &queue, oneapi::mkl::uplo upper_lower, std::int64_t n, sycl::buffer<T,1> &a, std::int64_t lda, sycl::buffer<realT,1> &d, sycl::buffer<realT,1> &e, sycl::buffer<T,1> &tau, sycl::buffer<T,1> &scratchpad, std::int64_t scratchpad_size)
}
Input Parameters
- queue
The queue where the routine should be executed.
- upper_lower
Must be
uplo::upperoruplo::lower.If
upper_lower = uplo::upper,astores the upper triangular part of \(A\).If
upper_lower = uplo::lower,astores the lower triangular part of \(A\).- n
The order of the matrices \(A\) \((0 \le n)\).
- a
Buffer, size
(lda,*). The bufferacontains either the upper or lower triangle of the Hermitian matrix \(A\), as specified by upper_lower.The second dimension of
amust be at least \(\max(1, n)\).- lda
The leading dimension of
a; at least \(\max(1, n)\)- scratchpad_size
Size of scratchpad memory as a number of floating point elements of type
T. Size should not be less than the value returned by hetrd_scratchpad_size function.
Output Parameters
- a
On exit,
if
upper_lower = uplo::upper, the diagonal and first superdiagonal of \(A\) are overwritten by the corresponding elements of the tridiagonal matrix \(T\), and the elements above the first superdiagonal, with the buffertau, represent the orthogonal matrix \(Q\) as a product of elementary reflectors;if
upper_lower = uplo::lower, the diagonal and first subdiagonal of \(A\) are overwritten by the corresponding elements of the tridiagonal matrix \(T\), and the elements below the first subdiagonal, with the buffertau, represent the orthogonal matrix \(Q\) as a product of elementary reflectors.- d
Buffer containing the diagonal elements of the matrix \(T\). The dimension of
dmust be at least \(\max(1, n)\).- e
Buffer containing the off diagonal elements of the matrix \(T\). The dimension of
emust be at least \(\max(1, n-1)\).- tau
Buffer, size at least \(\max(1, n-1)\). Stores \((n-1)\) scalars that define elementary reflectors in decomposition of the unitary matrix \(Q\) in a product of \(n-1\) elementary reflectors.
- scratchpad
Buffer holding scratchpad memory to be used by routine for storing intermediate results.
hetrd (USM Version)#
Syntax
namespace oneapi::mkl::lapack {
sycl::event hetrd(sycl::queue &queue, oneapi::mkl::uplo upper_lower, std::int64_t n, T *a, std::int64_t lda, RealT *d, RealT *e, T *tau, T *scratchpad, std::int64_t scratchpad_size, const std::vector<sycl::event> &events = {})
}
Input Parameters
- queue
The queue where the routine should be executed.
- upper_lower
Must be
uplo::upperoruplo::lower.If
upper_lower = uplo::upper,astores the upper triangular part of \(A\).If
upper_lower = uplo::lower,astores the lower triangular part of \(A\).- n
The order of the matrices \(A\) \((0 \le n)\).
- a
The pointer to matrix \(A\), size
(lda,*). Contains either the upper or lower triangle of the Hermitian matrix \(A\), as specified byupper_lower. The second dimension ofamust be at least \(\max(1, n)\).- lda
The leading dimension of
a; at least \(\max(1, n)\)- scratchpad_size
Size of scratchpad memory as a number of floating point elements of type
T. Size should not be less than the value returned by hetrd_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- a
On exit,
if
upper_lower = uplo::upper, the diagonal and first superdiagonal of \(A\) are overwritten by the corresponding elements of the tridiagonal matrix \(T\), and the elements above the first superdiagonal, with the arraytau, represent the orthogonal matrix \(Q\) as a product of elementary reflectors;if
upper_lower = uplo::lower, the diagonal and first subdiagonal of \(A\) are overwritten by the corresponding elements of the tridiagonal matrix \(T\), and the elements below the first subdiagonal, with the arraytau, represent the orthogonal matrix \(Q\) as a product of elementary reflectors.- d
Pointer to diagonal elements of the matrix \(T\). The dimension of
dmust be at least \(\max(1, n)\).- e
Pointer to off diagonal elements of the matrix \(T\). The dimension of
emust be at least \(\max(1, n-1)\).- tau
Pointer to array of size at least \(\max(1, n-1)\). Stores \((n-1)\) scalars that define elementary reflectors in decomposition of the unitary matrix \(Q\) in a product of \(n-1\) elementary reflectors.
- scratchpad
Pointer to scratchpad memory to be used by routine for storing intermediate results.
Return Values
Output event to wait on to ensure computation is complete.
Parent topic: LAPACK Singular Value and Eigenvalue Problem Routines