geqrf#
Computes the QR factorization of a general \(m \times n\) matrix.
Description
geqrf
supports the following precisions:
T |
---|
|
|
|
|
The routine forms the QR factorization of a general \(m \times n\) matrix \(A\). No pivoting is performed.
The routine does not form the matrix \(Q\) explicitly. Instead, \(Q\) is represented as a product of \(\min(m, n)\) elementary reflectors. Routines are provided to work with \(Q\) in this representation.
geqrf (Buffer Version)#
Syntax
namespace oneapi::mkl::lapack {
void geqrf(sycl::queue &queue, std::int64_t m, std::int64_t n, 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.
- m
The number of rows in the matrix \(A\) (\(0 \le m\)).
- n
The number of columns in \(A\) (\(0 \le n\)).
- a
Buffer holding input matrix \(A\). Must have size at least \(\text{lda} \cdot n\).
- lda
The leading dimension of \(A\); at least \(\max(1, m)\).
- 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 geqrf_scratchpad_size function.
Output Parameters
- a
Output buffer, overwritten by the factorization data as follows:
The elements on and above the diagonal of the array contain the \(\min(m,n) \times n\) upper trapezoidal matrix \(R\) (\(R\) is upper triangular if \(m \ge n\)); the elements below the diagonal, with the array tau, represent the orthogonal matrix \(Q\) as a product of \(\min(m,n)\) elementary reflectors.
- tau
Output buffer, size at least \(\max(1, \min(m, n))\). Contains scalars that define elementary reflectors for the matrix \(Q\) in its decomposition in a product of elementary reflectors.
- scratchpad
Buffer holding scratchpad memory to be used by routine for storing intermediate results.
geqrf (USM Version)#
Syntax
namespace oneapi::mkl::lapack {
sycl::event geqrf(sycl::queue &queue, std::int64_t m, std::int64_t n, 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.
- m
The number of rows in the matrix \(A\) (\(0 \le m\)).
- n
The number of columns in \(A\) (\(0 \le n\)).
- a
Pointer to memory holding input matrix \(A\). Must have size at least \(\text{lda} \cdot n\).
- lda
The leading dimension of \(A\); at least \(\max(1, m)\).
- 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 geqrf_scratchpad_size function.- events
List of events to wait for before starting computation. Defaults to empty list.
Output Parameters
- a
Overwritten by the factorization data as follows:
The elements on and above the diagonal of the array contain the \(\min(m,n) \times n\) upper trapezoidal matrix \(R\) (\(R\) is upper triangular if \(m \ge n\)); the elements below the diagonal, with the array tau, represent the orthogonal matrix \(Q\) as a product of \(\min(m,n)\) elementary reflectors.
- tau
Array, size at least \(\max(1, \min(m, n))\). Contains scalars that define elementary reflectors for the matrix \(Q\) in its decomposition in a product of elementary reflectors.
- 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 Linear Equation Routines