geqrf_batch#
Computes the QR factorizations of a batch of general matrices.
Description
geqrf_batch supports the following precisions.
T
float
double
std::complex<float>
std::complex<double>
geqrf_batch (Buffer Version)#
Description
The buffer version of geqrf_batch supports only the strided API.
Strided API
Syntax
namespace oneapi::mkl::lapack {
void geqrf_batch(sycl::queue &queue, std::int64_t m, std::int64_t n, sycl::buffer<T> &a, std::int64_t lda, std::int64_t stride_a, sycl::buffer<T> &tau, std::int64_t stride_tau, std::int64_t batch_size, sycl::buffer<T> &scratchpad, std::int64_t scratchpad_size)
}
Input Parameters
- queue
Device queue where calculations will be performed.
- m
Number of rows in matrices \(A_i\) (\(0 \le m\)).
- n
Number of columns in matrices \(A_i\) (\(0 \le n\)).
- a
Array holding input matrices \(A_i\).
- lda
Leading dimension of matrices \(A_i\).
- stride_a
Stride between the beginnings of matrices \(A_i\) inside the batch array
a.- stride_tau
Stride between the beginnings of arrays \(\tau_i\) inside the array
tau.- batch_size
Number of problems in a batch.
- scratchpad
Scratchpad memory to be used by routine for storing intermediate results.
- scratchpad_size
Size of scratchpad memory as the number of floating point elements of type
T. Size should not be less than the value returned by the Strided API of the geqrf_batch_scratchpad_size function.
Output Parameters
- a
Factorization data as follows: The elements on and above the diagonal of \(A_i\) contain the \(\min(m,n) \times n\) upper trapezoidal matrices \(R_i\) (\(R_i\) is upper triangular if \(m \ge n\)); the elements below the diagonal, with the array \(\tau_i\), contain the orthogonal matrix \(Q_i\) as a product of \(\min(m,n)\) elementary reflectors.
- tau
Array to store batch of \(\tau_i\), each of size \(\min(m,n)\), containing scalars that define elementary reflectors for the matrices \(Q_i\) in its decomposition in a product of elementary reflectors.
geqrf_batch (USM Version)#
Description
The USM version of geqrf_batch supports the group API and strided API.
Group API
The routine forms the \(Q_iR_i\) factorizations of a general \(m \times n\) matrices \(A_i\), \(i \in \{1...batch\_size\}\), where batch_size is the sum of all parameter group sizes as provided with group_sizes array.
No pivoting is performed during factorization.
The routine does not form the matrices \(Q_i\) explicitly. Instead, \(Q_i\) is represented as a product of \(\min(m,n)\) elementary reflectors. Routines are provided to work with \(Q_i\) in this representation.
The total number of problems to solve, batch_size, is a sum of sizes of all of the groups of parameters as provided by group_sizes array.
Syntax
namespace oneapi::mkl::lapack {
sycl::event geqrf_batch(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<sycl::event> &events = {})
}
Input Parameters
- queue
Device queue where calculations will be performed.
- m
Array of
group_count\(m_g\) parameters. Each \(m_g\) specifies the number of rows in matrices \(A_i\) from arraya, belonging to group \(g\).- n
Array of
group_count\(n_g\) parameters. Each \(n_g\) specifies the number of columns in matrices \(A_i\) from arraya, belonging to group \(g\).- a
Array of
batch_sizepointers to input matrices \(A_i\), each of size \(\text{lda}_g\cdot n_g\) (\(g\) is an index of group to which \(A_i\) belongs)- lda
Array of
group_count\(\text{lda}_g`\) parameters, each representing the leading dimensions of input matrices \(A_i\) from arraya, belonging to group \(g\).- group_count
Specifies the number of groups of parameters. Must be at least 0.
- group_sizes
Array of
group_countintegers. 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 the number of floating point elements of type
T. Size should not be less than the value returned by the Group API of the geqrf_batch_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- a
Factorization data as follows: The elements on and above the diagonal of \(A_i\) contain the \(\min(m_g,n_g) \times n_g\) upper trapezoidal matrices \(R_i\) (\(R_i\) is upper triangular if \(m_g \ge n_g\)); the elements below the diagonal, with the array \(\tau_i\), contain the orthogonal matrix \(Q_i\) as a product of \(\min(m_g,n_g)\) elementary reflectors. Here \(g\) is the index of the parameters group corresponding to the \(i\)-th decomposition.
- tau
Array of pointers to store arrays \(\tau_i\), each of size \(\min(m_g,n_g)\), containing scalars that define elementary reflectors for the matrices \(Q_i\) in its decomposition in a product of elementary reflectors. Here \(g\) is the index of the parameters group corresponding to the \(i\)-th decomposition.
Return Values
Output event to wait on to ensure computation is complete.
Strided API
The routine forms the \(Q_iR_i\) factorizations of general \(m \times n\) matrices \(A_i\). No pivoting is performed. The routine does not form the matrices \(Q_i\) explicitly. Instead, \(Q_i\) is represented as a product of \(\min(m,n)\) elementary reflectors. Routines are provided to work with \(Q_i\) in this representation.
Syntax
namespace oneapi::mkl::lapack {
sycl::event geqrf_batch(sycl::queue &queue, std::int64_t m, std::int64_t n, T *a, std::int64_t lda, std::int64_t stride_a, T *tau, std::int64_t stride_tau, std::int64_t batch_size, T *scratchpad, std::int64_t scratchpad_size, const std::vector<sycl::event> &events = {})
}
Input Parameters
- queue
Device queue where calculations will be performed.
- m
Number of rows in matrices \(A_i\) (\(0 \le m\)).
- n
Number of columns in matrices \(A_i\) (\(0 \le n\)).
- a
Array holding input matrices \(A_i\).
- lda
Leading dimensions of \(A_i\).
- stride_a
Stride between the beginnings of matrices \(A_i\) inside the batch array
a.- stride_tau
Stride between the beginnings of arrays \(\tau_i\) inside the array
tau.- batch_size
Number of problems in a batch.
- scratchpad
Scratchpad memory to be used by routine for storing intermediate results.
- scratchpad_size
Size of scratchpad memory as the number of floating point elements of type
T. Size should not be less than the value returned by the Strided API of the geqrf_batch_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- a
Factorization data as follows: The elements on and above the diagonal of \(A_i\) contain the \(\min(m,n) \times n\) upper trapezoidal matrices \(R_i\) (\(R_i\) is upper triangular if \(m \ge n\)); the elements below the diagonal, with the array \(\tau_i\), contain the orthogonal matrix \(Q_i\) as a product of \(\min(m,n)\) elementary reflectors.
- tau
Array to store batch of \(\tau_i\), each of size \(\min(m,n)\), containing scalars that define elementary reflectors for the matrices \(Q_i\) in its decomposition in a product of elementary reflectors.
Return Values
Output event to wait on to ensure computation is complete.
Parent topic: LAPACK-like Extensions Routines