ungtr¶
Generates the complex unitary matrix Q determined by
hetrd. This routine
belongs to the oneapi::mkl::lapack
namespace.
Syntax
namespace oneapi::mkl::lapack {
void ungtr(cl::sycl::queue &queue, mkl::uplo uplo, std::int64_t n, cl::sycl::buffer<T> &a, std::int64_t lda, cl::sycl::buffer<T> &tau, cl::sycl::buffer<T> &scratchpad, std::int64_t scratchpad_size)
}
ungtr
supports the following precision and devices.
T |
Devices Supported |
---|---|
|
Host and CPU |
|
Host and CPU |
Description
The routine explicitly generates the n
-by-n
unitary matrix
Q
formed by
hetrd when
reducing a complex Hermitian matrix A
to tridiagonal form:
A = Q*T*QH
. Use this routine after a call to
hetrd.
Input Parameters
- queue
Device queue where calculations will be performed.
- uplo
Must be uplo::upper or uplo::lower. Uses the same uplo as supplied to the hetrd function.
- n
The order of matrix
Q
(0≤n
).- a
The buffer
a
returned by the hetrd function. The second dimension ofa
must be at leastmax(1, n)
.- lda
The leading dimension of
a
(n≤lda
).- tau
The buffer tau returned by the hetrd function. The dimension of tau must be at least
max(1,n-1)
.- scratchpad
Buffer holding scratchpad memory to be used by the routine for storing intermediate results.
- 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 the ungtr_scratchpad_size function.
Output Parameters
- a
Overwritten by the unitary matrix
Q
.
Exceptions
mkl::lapack::exception |
This exception is thrown when problems occur during calculations. You can obtain the info code of the problem using the info() method of the exception object: If |