ungqr_batch (Group Version)#
Generates the complex unitary matrices Qi of the batch of QR
factorizations formed by geqrf_batch. This routine belongs to the
oneapi::mkl::lapack namespace.
Description#
The routine generates the wholes or parts of m-by-m unitary
matrix Qi of the QR factorization formed by the
geqrf_batch (Group Version)
function.
Usually Qi is determined from the QR factorization
of an m by p matrix Ai with m ≥ p. To compute
the whole matrices Qi, use:
ungqr_batch(queue, m, m, p, a, ...)
To compute the leading p columns of Qi (which form
an orthonormal basis in the space spanned by the columns of
Ai):
ungqr_batch(queue, m, p, p, a, ...)
To compute the matrix Qik of the QR
factorization of the leading k columns of the matrices
Ai:
ungqr_batch(queue, m, m, k, a, ...)
To compute the leading k columns of Qik
(which form an orthonormal basis in the space spanned by the leading
k columns of the matrices Ai):
ungqr_batch(queue, m, k, k, a, ...)
API#
Syntax#
namespace oneapi::mkl::lapack {
sycl::event ungqr_batch(sycl::queue &queue,
int64_t *m,
int64_t *n,
int64_t *k,
T **a,
int64_t *lda,
const T * const *tau,
int64_t group_count,
int64_t *group_sizes,
T *scratchpad,
int64_t scratchpad_size,
const std::vector<sycl::event> &events = {})
}
This function supports the following precisions and devices:
T |
Devices supported |
|---|---|
|
CPU and GPU* |
|
CPU and GPU* |
*Interface support only; all computations are performed on the CPU.
Input Parameters#
- queue
Device queue where calculations will be performed.
- m
Array of
group_countparameters mg as previously supplied to geqrf_batch (Group Version) (mg ≥ 0).- n
Array of
group_countparameters ng as previously supplied to geqrf_batch (Group Version) (ng ≥ 0).- k
Array of
group_countparameterskg. The number of elementary reflectors whose product defines the matricesQi (0 ≤kg ≤ng) .- a
Array resulting after a call to geqrf_batch (Group Version).
- lda
Array of leading dimension of
Ai as previously supplied to geqrf_batch (Group Version) (ldag ≥max(1, mg)).- tau
Array resulting after a call to geqrf_batch (Group Version).
- group_count
Specifies the number of groups of parameters. Must be at least 0.
- group_sizes
Array of group_count integers. Array element with index
gspecifies the number of problems to solve for each of the groups of parametersg. So the total number of problems to solve,batch_size, is a sum of all parameter group sizes.- scratchpad
Scratchpad memory to be used by 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 then the value returned by ungqr_batch_scratchpad_size (Group Version).
- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters#
- a
Matrices pointed to by array
aare overwritten byng leading columns of themg-by-mg unitary matricesQi, wheregis an index of group of parameters corresponding toQi.
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.