ungbr (USM Version)¶
Generates the complex unitary matrix Q or Pt determined by the
gebrd
(USM Version) function. This routine belongs to the oneapi::mkl::lapack
namespace.
Description¶
The routine generates the whole or part of the unitary matrices Q
and P
T formed by the gebrd (USM Version)
function. All valid combinations of arguments are described in Input
Parameters; in most cases you need the following:
To compute the whole m
-by-m
matrix Q
, use:
ungbr(queue, generate::q, m, m, n, a, ...)
(note that the array a
must have at least m
columns).
To form the n
leading columns of Q
if m > n
, use:
ungbr(queue, generate::q, m, n, n, a, ...)
To compute the whole n
-by-n
matrix P
T, use:
ungbr(queue, generate::p, n, n, m, a, ...)
(note that the array a
must have at least n
rows).
To form the m
leading rows of P
T if m < n
, use:
ungbr(queue, generate::p, m, n, m, a, ...)
API¶
Syntax¶
namespace oneapi::mkl::lapack {
cl::sycl::event ungbr(cl::sycl::queue &queue,
mkl::generate gen,
std::int64_t m,
std::int64_t n,
std::int64_t k,
T *a,
std::int64_t lda,
T *tau,
T *scratchpad,
std::int64_t scratchpad_size,
const std::vector<cl::sycl::event> &events = {})
}
ungbr
(USM version) supports the following precisions and
devices:
T |
Devices supported |
---|---|
|
Host and CPU |
|
Host and CPU |
Input Parameters¶
- queue
Device queue where calculations will be performed.
- gen
Must be
generate::q
orgenerate::p
.If gen
= generate::q
, the routine generates the matrixQ
.If gen
= generate::p
, the routine generates the matrixP
T.- m
The number of rows in the matrix
Q
orP
T to be returned(0≤m)
.If gen
= generate::q
,m ≥ n ≥ min(m, k)
.If gen
= generate::p
,n ≥ m ≥ min(n, k)
.- n
The number of rows in the matrix
Q
orP
T to be returned(0≤n)
. See m for constraints.- k
If gen
= generate::q
, the number of columns in the originalm
-by-k matrix reduced by gebrd (USM Version).If gen
= generate::p
, the number of rows in the originalk
-by-n matrix reduced by gebrd (USM Version).- a
Pointer to memory returned by gebrd (USM Version).
- lda
The leading dimension of a.
- tau
For gen
= generate::q
, the array tauq is returned by the gebrd (USM Version) function.For gen
= generate::p
, the array taup is returned by the gebrd (USM Version) function.The dimension of tau must be at least
min(m,k)
forgen = generate::q
, ormin(n,k)
forgen = generate::p
.- 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 ungbr_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters¶
- a
Overwritten by n leading columns of the m-by-m orthogonal matrix
Q
orP
T (or the leading rows or columns thereof) as specified by gen, m, and n.
Exceptions¶
Exception |
Description |
---|---|
|
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 If |
Return Values¶
Output event to wait on to ensure computation is complete.