unmtr (USM Version)¶
Multiplies a complex matrix by the complex unitary matrix Q determined
by the hetrd (USM Version) function. This routine belongs to the
oneapi::mkl::lapack
namespace.
Syntax
-
cl::sycl::event
unmtr
(cl::sycl::queue &queue, mkl::side side, mkl::uplo uplo, mkl::transpose trans, std::int64_t m, std::int64_t n, T *a, std::int64_t lda, T *tau, T *c, std::int64_t ldc, T *scratchpad, std::int64_t scratchpad_size, const cl::sycl::vector_class<cl::sycl::event> &events = {})¶
unmtr
(USM version) supports the following precisions and
devices:
T |
Devices supported |
---|---|
|
Host and CPU |
|
Host and CPU |
Description
The routine multiplies a complex matrix C
by Q
or
Q
H, where Q
is the unitary matrix Q
formed by
hetrd (USM
Version) when
reducing a complex Hermitian matrix A
to tridiagonal form:
A = Q*T*QT
. Use this routine after a call to hetrd (USM
Version).
Depending on the parameters side and trans, the routine can form one
of the matrix products Q*C
, Q
T*C
, C*Q
, or
C*Q
T (overwriting the result on C
).
Input Parameters
In the descriptions below, r
denotes the order of Q
:
|
if |
---|---|
|
if |
- queue
Device queue where calculations will be performed.
- side
Must be either
side::left
orside::right
.If
side = side::left
,Q
orQ
T is applied toC
from the left.If
side = side::right
,Q
orQ
T is applied toC
from the right.- uplo
Must be either
uplo::upper
oruplo::lower
. Uses the sameuplo
as supplied to hetrd (USM Version).- trans
Must be either
transpose::nontrans
ortranspose::trans
.If
trans = transpose::nontrans
, the routine multipliesC
byQ
.If
trans = transpose::trans
, the routine multipliesC
byQ
T.- m
The number of rows in the matrix
C
(m ≥ 0)
.- n
The number of columns in the matrix
C
(n ≥ 0)
.- a
The pointer to memory returned by hetrd (USM Version).
- lda
The leading dimension of a
(max(1, r) ≤ lda)
.- tau
The pointer to tau returned by hetrd (USM Version). The dimension of tau must be at least
max(1, r-1)
.- c
The pointer to memory containing the matrix
C
. The second dimension of c must be at leastmax(1, n)
.- ldc
The leading dimension of c
(max(1, n) ≤ ldc)
.- scratchpad
Pointer to 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 unmtr_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- c
Overwritten by the product
Q*C
,Q
T*C
,C*Q
, orC*Q
T (as specified by side and trans).
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 |
Return Values
Output event to wait on to ensure computation is complete.