ungbr#
Generates the complex unitary matrix \(Q\) or \(P^{t}\) determined by gebrd.
Description
ungbr
supports the following precisions.
T
std::complex<float>
std::complex<double>
The routine generates the whole or part of the unitary matrices \(Q\) and \(P^{H}\) 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 \times m\) matrix \(Q\), use:
oneapi::mkl::lapack::ungbr(queue, generate::q, m, m, n, a, ...)
(note that the buffer a
must have at least \(m\) columns).
To form the \(n\) leading columns of \(Q\) if \(m > n\), use:
oneapi::mkl::lapack::ungbr(queue, generate::q, m, n, n, a, ...)
To compute the whole \(n \times n\) matrix \(P^{T}\), use:
oneapi::mkl::lapack::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:
oneapi::mkl::lapack::ungbr(queue, generate::p, m, n, m, a, ...)
ungbr (Buffer Version)#
Syntax
namespace oneapi::mkl::lapack {
void ungbr(sycl::queue &queue, oneapi::mkl::generate gen, std::int64_t m, std::int64_t n, std::int64_t k, sycl::buffer<T,1> &a, std::int64_t lda, sycl::buffer<T,1> &tau, 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::q
orgenerate::p
.If
gen = generate::q
, the routine generates the matrix \(Q\).If
gen = generate::p
, the routine generates the matrix \(P^{T}\).- m
The number of rows in the matrix \(Q\) or \(P^{T}\) to be returned \((0 \le m)\).
If
gen = generate::q
, \(m \ge n \ge \min(m, k)\).If
gen = generate::p
, \(n \ge m \ge \min(n, k)\).- n
The number of columns in the matrix \(Q\) or \(P^{T}\) to be returned \((0 \le n)\). See
m
for constraints.- k
If
gen = generate::q
, the number of columns in the original \(m \times k\) matrix returned by gebrd.If
gen = generate::p
, the number of rows in the original \(k \times n\) matrix returned by gebrd.- a
The buffer
a
as returned by gebrd.- lda
The leading dimension of
a
.- tau
For
gen = generate::q
, the arraytauq
as returned by gebrd. Forgen = generate::p
, the arraytaup
as returned by gebrd.The dimension of
tau
must be at least \(\max(1, \min(m, k))\) forgen = generate::q
, or \(\max(1, \min(m, k))\) forgen = 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 \times m\) unitary matrix \(Q\) or \(P^{T}\), (or the leading rows or columns thereof) as specified by
gen
,m
, andn
.- scratchpad
Buffer holding scratchpad memory to be used by routine for storing intermediate results.
ungbr (USM Version)#
Syntax
namespace oneapi::mkl::lapack {
sycl::event ungbr(sycl::queue &queue, oneapi::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<sycl::event> &events = {})
}
Input Parameters
- queue
The queue where the routine should be executed.
- gen
Must be
generate::q
orgenerate::p
.If
gen = generate::q
, the routine generates the matrix \(Q\).If
gen = generate::p
, the routine generates the matrix \(P^{T}\).- m
The number of rows in the matrix \(Q\) or \(P^{T}\) to be returned \((0 \ge m)\).
If
gen = generate::q
, \(m \ge n \ge \min(m, k)\).If
gen = generate::p
, \(n \ge m \ge \min(n, k)\).- n
The number of columns in the matrix \(Q\) or \(P^{T}\) to be returned \((0 \le n)\). See
m
for constraints.- k
If
gen = generate::q
, the number of columns in the original \(m \times k\) matrix returned by gebrd.If
gen = generate::p
, the number of rows in the original \(k \times n\) matrix returned by gebrd.- a
The pointer to
a
as returned by gebrd.- lda
The leading dimension of
a
.- tau
For
gen = generate::q
, the arraytauq
as returned by gebrd. Forgen = generate::p
, the arraytaup
as returned by gebrd.The dimension of
tau
must be at least \(\max(1, \min(m, k))\) forgen = generate::q
, or \(\max(1, \min(m, k))\) forgen = 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 \times m\) unitary matrix \(Q\) or \(P^{T}\), (or the leading rows or columns thereof) as specified by
gen
,m
, andn
.- scratchpad
Pointer to scratchpad memory to be used by routine for storing intermediate results.
Return Values
Output event to wait on to ensure computation is complete.
Parent topic: LAPACK Singular Value and Eigenvalue Problem Routines