geqrf_batch (Group Version)¶
Computes the batch of QR factorizations of a general m-by-n matrices.
This routine belongs to the oneapi::mkl::lapack namespace.
Description¶
The routine forms the QiRi factorizations
of a general m-by-n matrix Ai,
iϵ{1...batch_size}, where batch_size is a sum of all
parameter group sizes as provided with the group_sizes array. No
pivoting is performed during factorization.
The routine does not form the matrix Qi explicitly.
Instead, Qi is represented as a product of min(m,
n) elementary reflectors. Routines are provided to work with
Qi in this representation.
Total number of problems to solve, batch_size, is a sum of sizes
of all of the groups of parameters as provided by
group_sizesarray.
API¶
Syntax¶
namespace oneapi::mkl::lapack {
cl::sycl::event geqrf_batch(cl::sycl::queue &queue,
std::int64_t *m,
std::int64_t *n,
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<cl::sycl::event> &events = {})
}
Function supports the following precisions and devices.
T
Devices supported
floatHost, CPU, and GPU
doubleHost, CPU, and GPU
std::complex<float>Host, CPU, and GPU
std::complex<double>Host, CPU, and GPU
Input Parameters¶
- queue
Device queue where calculations will be performed.
- m
Array of
group_countparameters mg parameters.Each of mg specifies the number of rows in the matrices
Ai from arraya, belonging to groupg.- n
Array of
group_countparameters ng parameters.Each of ng specifies the number of columns in the matrices
Ai from arraya, belonging to groupg.- a
Array of
batch_sizepointers to input matricesAi, each being of sizeldag*ng (gis an index of group to whichAi belongs).- lda
Array of
group_countldag parameters, each representing the leading dimensions of input matricesAi, from arraya, belonging to groupg.- 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 geqrf_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 by the factorization data as follows:The elements on and above the diagonal of
Ai contain themin(mg,mng)-by-ng upper trapezoidal matricesRi (Ri is upper triangular ifmg≥ng); the elements below the diagonal, with the arraytaui, present the orthogonal matrixQi as a rpoduct ofmin(mg,ng)elementary reflectors.Here,
gis an index of parameters group corresponding toi-th decomposition.- tau
Array of pointers to store
taui, each of sizemg,ng, containing scalars that define elementary reflectors for the matricesQi in its decomposition in a product of elementary reflectors.Here,
gis an index of parameters group corresponding toi-th decomposition.
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.