ungbr¶
Generates the complex unitary matrix Q or Pt determined by gebrd.
ungbrsupports the following precisions.
T
std::complex<float>
std::complex<double>
Description
The routine generates the whole or part of the unitary matrices Q
and PH formed by the routines
gebrd.
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 buffera 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 PT, 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 PT if m < n, use:
ungbr(queue, generate::p, m, n, m, a, …)
ungbr (BUFFER Version)¶
Syntax
-
void
onemkl::lapack::ungbr(cl::sycl::queue &queue, onemkl::generate gen, std::int64_t m, std::int64_t n, std::int64_t k, cl::sycl::buffer<T, 1> &a, std::int64_t lda, cl::sycl::buffer<T, 1> &tau, cl::sycl::buffer<T, 1> &scratchpad, std::int64_t scratchpad_size)¶
Input Parameters
- queue
The queue where the routine should be executed.
- gen
Must be
generate::qorgenerate::p.If
gen = generate::q, the routine generates the matrixQ.If
gen = generate::p, the routine generates the matrixPT.- m
The number of rows in the matrix
QorPT 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 columns in the matrix
QorPT to be returned(0 ≤ n). See m for constraints.- k
If
gen = generate::q, the number of columns in the originalm-by-kmatrix returned by gebrd.If
gen = generate::p, the number of rows in the originalk-by-nmatrix returned by gebrd.- a
The buffer
aas returned by gebrd.- lda
The leading dimension of a.
- tau
For gen
= generate::q, the arraytauqas returned by gebrd. For gen= generate::p, the arraytaupas returned by gebrd.The dimension of tau must be at least
max(1, min(m, k))for gen=generate::q, ormax(1, min(m, k))for gen= generate::p.- 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 ungbr_scratchpad_size function.
Output Parameters
- a
Overwritten by n leading columns of the m-by-m unitary matrix
QorPT, (or the leading rows or columns thereof) as specified by gen, m, and n.- scratchpad
Buffer holding scratchpad memory to be used by routine for storing intermediate results.
Throws
- onemkl::lapack::exception
Exception is thrown in case of problems happened during calculations. The
infocode of the problem can be obtained by get_info() method of exception object:If
info=-i, thei-th parameter had an illegal value.If
infoequals to value passed as scratchpad size, andget_detail()returns non zero, then passed scratchpad is of insufficient size, and required size should not be less than value return byget_detail()method of exception object.
ungbr (USM Version)¶
Syntax
-
cl::sycl::event
onemkl::lapack::ungbr(cl::sycl::queue &queue, onemkl::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 cl::sycl::vector_class<cl::sycl::event> &events = {})¶
Input Parameters
- queue
The queue where the routine should be executed.
- gen
Must be
generate::qorgenerate::p.If
gen = generate::q, the routine generates the matrixQ.If
gen = generate::p, the routine generates the matrixPT.- m
The number of rows in the matrix
QorPT 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 columns in the matrix
QorPT to be returned(0 ≤ n). See m for constraints.- k
If
gen = generate::q, the number of columns in the originalm-by-kmatrix returned by gebrd.If
gen = generate::p, the number of rows in the originalk-by-nmatrix returned by gebrd.- a
The pointer to
aas returned by gebrd.- lda
The leading dimension of a.
- tau
For gen
= generate::q, the arraytauqas returned by gebrd. For gen= generate::p, the arraytaupas returned by gebrd.The dimension of tau must be at least
max(1, min(m, k))for gen=generate::q, ormax(1, min(m, k))for gen= generate::p.- 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 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 unitary matrix
QorPT, (or the leading rows or columns thereof) as specified by gen, m, and n.- scratchpad
Pointer to scratchpad memory to be used by routine for storing intermediate results.
Throws
- onemkl::lapack::exception
Exception is thrown in case of problems happened during calculations. The
infocode of the problem can be obtained by get_info() method of exception object:If
info=-i, thei-th parameter had an illegal value.If
infoequals to value passed as scratchpad size, andget_detail()returns non zero, then passed scratchpad is of insufficient size, and required size should not be less than value return byget_detail()method of exception object.
Return Values
Output event to wait on to ensure computation is complete.
Parent topic: LAPACK Singular Value and Eigenvalue Problem Routines