ungtr#
Generates the complex unitary matrix \(Q\) determined by hetrd.
Description
ungtr
supports the following precisions.
T
std::complex<float>
std::complex<double>
The routine explicitly generates the \(n \times n\) unitary matrix \(Q\) formed by hetrd when reducing a complex Hermitian matrix \(A\) to tridiagonal form: \(A = QTQ^H\). Use this routine after a call to hetrd.
ungtr (Buffer Version)#
Syntax
namespace oneapi::mkl::lapack {
void ungtr(sycl::queue &queue, oneapi::mkl::uplo upper_lower, std::int64_t n, sycl::buffer<T,1> &a, std::int64_t lda, 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::upper
oruplo::lower
. Uses the sameupper_lower
as supplied to hetrd.- n
The order of the matrix \(Q\) \((0 \le n)\).
- a
The buffer
a
as returned by hetrd. The second dimension ofa
must be at least \(\max(1, n)\).- lda
The leading dimension of
a
\((n \le \text{lda})\).- tau
The buffer
tau
as returned by hetrd. The dimension oftau
must be at least \(\max(1, n-1)\).- 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 ungtr_scratchpad_size function.
Output Parameters
- a
Overwritten by the unitary matrix \(Q\).
- scratchpad
Buffer holding scratchpad memory to be used by routine for storing intermediate results.
ungtr (USM Version)#
Syntax
namespace oneapi::mkl::lapack {
sycl::event ungtr(sycl::queue &queue, oneapi::mkl::uplo upper_lower, std::int64_t n, T *a, std::int64_t lda, 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::upper
oruplo::lower
. Uses the sameupper_lower
as supplied to hetrd.- n
The order of the matrix \(Q\) \((0 \le n)\).
- a
The pointer to
a
as returned by hetrd. The second dimension ofa
must be at least \(\max(1, n)\).- lda
The leading dimension of
a
\((n \le \text{lda})\).- tau
The pointer to
tau
as returned by hetrd. The dimension oftau
must be at least \(\max(1, n-1)\).- 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 ungtr_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- a
Overwritten by the unitary matrix \(Q\).
- 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