hetrd¶
Reduces a complex Hermitian matrix to tridiagonal form.
hetrd
supports the following precisions.
Routine name
T
chetrd
std::complex<float>
zhetrd
std::complex<double>
Description
The routine reduces a complex Hermitian matrix A
to symmetric
tridiagonal form T
by a unitary similarity transformation:
A = Q*T*QH
. 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
-
void
onemkl::lapack
::
hetrd
(cl::sycl::queue &queue, onemkl::uplo upper_lower, std::int64_t n, cl::sycl::buffer<T, 1> &a, std::int64_t lda, cl::sycl::buffer<realT, 1> &d, cl::sycl::buffer<realT, 1> &e, cl::sycl::buffer<T, 1> &tau, cl::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::upper
oruplo::lower
.If
upper_lower = uplo::upper
, a stores the upper triangular part ofA
.If
upper_lower = uplo::lower
, a stores the lower triangular part ofA
.- n
The order of the matrices
A
(0≤n)
.- a
Buffer, size
(lda,*)
. The buffer a contains either the upper or lower triangle of the Hermitian matrixA
, as specified by upper_lower.The second dimension of a must 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 ofA
are overwritten by the corresponding elements of the tridiagonal matrixT
, and the elements above the first superdiagonal, with the buffer tau, represent the orthogonal matrixQ
as a product of elementary reflectors;if
upper_lower = uplo::lower
, the diagonal and first subdiagonal ofA
are overwritten by the corresponding elements of the tridiagonal matrixT
, and the elements below the first subdiagonal, with the buffer tau, represent the orthogonal matrixQ
as a product of elementary reflectors.- d
Buffer containing the diagonal elements of the matrix
T
. The dimension of d must be at leastmax(1, n)
.- e
Buffer containing the off diagonal elements of the matrix
T
. The dimension of e must be at leastmax(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 matrixQ
in a product ofn-1
elementary reflectors.- scratchpad
Buffer holding scratchpad memory to be used by routine for storing intermediate results.
Throws
- onemkl::lapack::exception
Exception is thrown in case of problems happened during calculations. The
info
code of the problem can be obtained by get_info() method of exception object:If
info=-i
, thei
-th parameter had an illegal value.If
info
equals to value passed as scratchpad size, andget_detail()
returns non zero, then passed scratchpad is of insufficient size, and required size should not be less than value return byget_detail()
method of exception object.
hetrd (USM Version)¶
Syntax
-
cl::sycl::event
onemkl::lapack
::
hetrd
(cl::sycl::queue &queue, onemkl::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 cl::sycl::vector_class<cl::sycl::event> &events = {})¶
Input Parameters
- queue
The queue where the routine should be executed.
- upper_lower
Must be
uplo::upper
oruplo::lower
.If
upper_lower = uplo::upper
, a stores the upper triangular part ofA
.If
upper_lower = uplo::lower
, a stores the lower triangular part ofA
.- n
The order of the matrices
A
(0≤n)
.- a
The pointer to matrix
A
, size(lda,*)
. Contains either the upper or lower triangle of the Hermitian matrixA
, as specified byupper_lower
. The second dimension of a must be at leastmax(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 ofA
are overwritten by the corresponding elements of the tridiagonal matrixT
, and the elements above the first superdiagonal, with the array tau, represent the orthogonal matrixQ
as a product of elementary reflectors;if
upper_lower = uplo::lower
, the diagonal and first subdiagonal ofA
are overwritten by the corresponding elements of the tridiagonal matrixT
, and the elements below the first subdiagonal, with the array tau, represent the orthogonal matrixQ
as a product of elementary reflectors.- d
Pointer to diagonal elements of the matrix
T
. The dimension of d must be at leastmax(1, n)
.- e
Pointer to off diagonal elements of the matrix
T
. The dimension of e must be at leastmax(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 matrixQ
in a product ofn-1
elementary reflectors.- scratchpad
Pointer to scratchpad memory to be used by routine for storing intermediate results.
Throws
- onemkl::lapack::exception
Exception is thrown in case of problems happened during calculations. The
info
code of the problem can be obtained by get_info() method of exception object:If
info=-i
, thei
-th parameter had an illegal value.If
info
equals to value passed as scratchpad size, andget_detail()
returns non zero, then passed scratchpad is of insufficient size, and required size should not be less than value return byget_detail()
method of exception object.
Return Values
Output event to wait on to ensure computation is complete.
Parent topic: LAPACK Singular Value and Eigenvalue Problem Routines