Visible to Intel only — GUID: GUID-CEC25B18-76F1-43E7-8269-48F5D3992810
Visible to Intel only — GUID: GUID-CEC25B18-76F1-43E7-8269-48F5D3992810
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,
std::int64_t *m,
std::int64_t *n,
std::int64_t *k,
T **a,
std::int64_t *lda,
T **tau,
std::int64_t group_count,
std::int64_t *group_sizes,
T *scratchpad,
std::int64_t scratchpad_size,
const std::vector<sycl::event> &events = {})
}
This function supports the following precisions and devices:
T |
Devices supported |
---|---|
std::complex<float> |
CPU and GPU* |
std::complex<double> |
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_count parameters mg as previously supplied to geqrf_batch (Group Version) (mg ≥ 0).
- n
-
Array of group_count parameters ng as previously supplied to geqrf_batch (Group Version) (ng ≥ 0).
- k
-
Array of group_count parameters kg. The number of elementary reflectors whose product defines the matrices Qi (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 g specifies the number of problems to solve for each of the groups of parameters g. 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 a are overwritten by ng leading columns of the mg-by-mg unitary matrices Qi, where g is an index of group of parameters corresponding to Qi.
Exceptions
Exception |
Description |
---|---|
mkl::lapack::batch_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 info = -n, the n-th parameter had an illegal value. If info equals the value passed as scratchpad size, and detail() returns non-zero, then the passed scratchpad is of insufficient size, and the required size should be not less then value returned by the detail() method of the exception object. |
Return Values
Output event to wait on to ensure computation is complete.