LQ Factorization#
GELQF (GEneral LQ Factorization) function computes batched LQ factorization of a general matrix:
using Householder reflection transformations, where:
A
is a batchedM x N
matrix, with leading dimensionlda >= M
if matrixA
is in column-major layout, orlda >= N
if matrixA
is row-major,Q
is an unitaryN x N
matrix, andL
is a lower triangular matrix ifM <= N
, or lower trapezoidal matrix ifM > N
.
The LQ factorization of A
is essentially the same as the QR factorization of \(A^T\), or \(A^H\) for complex data type.
cuSolverDx gelqf
device functions are (see Execution Methods):
__device__ void execute(data_type* A, data_type* tau);
// with runtime lda
__device__ void execute(data_type* A, const unsigned int lda, data_type* tau);
After the function, the lower triangular or lower trapezoidal part of input A
, including diagonal elements, is replaced by the matrix L
.
Matrix Q
is not explicitly formed. The elements above the diagonal of A
, with the array tau
, represent the matrix Q
as a product of min(M, N)
Householder vectors:
Each Householder vector has the form \(H(i) = I - tau[i] * v * v^H\), where:
tau
is an array of sizemin(M, N)
for each batch, andv
is a vector of sizeN
for each batch, withv[0 : i - 1] = 0
,v[i] = 1
, andconjugate(v[i+1 : N])
is stored on exit inA[i, i+1 : N]
.
The functions support A
being either column- or row-major memory layout, see Arrangement operator.