Generates the complex unitary matrix Q determined by hetrd. This routine belongs to the mkl::lapack namespace.
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 |
---|---|
std::complex<float> | Host and CPU |
std::complex<double> | Host and CPU |
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.
Device queue where calculations will be performed.
Must be uplo::upper or uplo::lower. Uses the same uplo as supplied to the hetrd function.
The order of matrix Q (0≤n).
The buffer a returned by the hetrd function. The second dimension of a must be at least max(1, n).
The leading dimension of a (n≤lda).
The buffer tau returned by the hetrd function. The dimension of tau must be at least max(1,n-1).
Buffer holding scratchpad memory to be used by the routine for storing intermediate results.
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.
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 info = -i, the i-th parameter had an illegal value. If info is equal to the value passed as scratchpad size, and detail() returns non-zero, then the passed scratchpad has an insufficient size, and the required size should not be less than the value returned by the detail() method of the exception object. |