Developer Reference for Intel® oneAPI Math Kernel Library for C

ID 766684
Date 12/16/2022
Public

A newer version of this document is available. Customers should click here to go to the newest version.

Document Table of Contents

?gels

Uses QR or LQ factorization to solve a overdetermined or underdetermined linear system with full rank matrix.

Syntax

lapack_int LAPACKE_sgels (int matrix_layout, char trans, lapack_int m, lapack_int n, lapack_int nrhs, float* a, lapack_int lda, float* b, lapack_int ldb);

lapack_int LAPACKE_dgels (int matrix_layout, char trans, lapack_int m, lapack_int n, lapack_int nrhs, double* a, lapack_int lda, double* b, lapack_int ldb);

lapack_int LAPACKE_cgels (int matrix_layout, char trans, lapack_int m, lapack_int n, lapack_int nrhs, lapack_complex_float* a, lapack_int lda, lapack_complex_float* b, lapack_int ldb);

lapack_int LAPACKE_zgels (int matrix_layout, char trans, lapack_int m, lapack_int n, lapack_int nrhs, lapack_complex_double* a, lapack_int lda, lapack_complex_double* b, lapack_int ldb);

Include Files
  • mkl.h
Description

The routine solves overdetermined or underdetermined real/ complex linear systems involving an m-by-n matrix A, or its transpose/ conjugate-transpose, using a QR or LQ factorization of A. It is assumed that A has full rank.

The following options are provided:

1. If trans = 'N' and mn: find the least squares solution of an overdetermined system, that is, solve the least squares problem

minimize ||b - A*x||2

2. If trans = 'N' and m < n: find the minimum norm solution of an underdetermined system A*X = B.

3. If trans = 'T' or 'C' and mn: find the minimum norm solution of an undetermined system AH*X = B.

4. If trans = 'T' or 'C' and m < n: find the least squares solution of an overdetermined system, that is, solve the least squares problem

minimize ||b - AH*x||2

Several right hand side vectors b and solution vectors x can be handled in a single call; they are formed by the columns of the right hand side matrix B and the solution matrix X (when coefficient matrix is A, B is m-by-nrhs and X is n-by-nrhs; if the coefficient matrix is AT or AH, B isn-by-nrhs and X is m-by-nrhs.

Input Parameters
matrix_layout

Specifies whether matrix storage layout is row major (LAPACK_ROW_MAJOR) or column major (LAPACK_COL_MAJOR).

trans

Must be 'N', 'T', or 'C'.

If trans = 'N', the linear system involves matrix A;

If trans = 'T', the linear system involves the transposed matrix AT (for real flavors only);

If trans = 'C', the linear system involves the conjugate-transposed matrix AH (for complex flavors only).

m

The number of rows of the matrix A (m 0).

n

The number of columns of the matrix A

(n 0).

nrhs

The number of right-hand sides; the number of columns in B (nrhs 0).

a, b

Arrays:

a(size max(1, lda*n) for column major layout and max(1, lda*m) for row major layout) contains the m-by-n matrix A.

b(size max(1, ldb*nrhs) for column major layout and max(1, ldb*max(m, n)) for row major layout) contains the matrix B of right hand side vectors.

lda

The leading dimension of a; at least max(1, m) for column major layout and at least max(1, n) for row major layout.

ldb

The leading dimension of b; must be at least max(1, m, n) for column major layout if trans='N' and at least max(1, n) if trans='T' and at least max(1, nrhs) for row major layout regardless of the value of trans.

Output Parameters
a

On exit, overwritten by the factorization data as follows:

if mn, array a contains the details of the QR factorization of the matrix A as returned by ?geqrf;

if m < n, array a contains the details of the LQ factorization of the matrix A as returned by ?gelqf.

b

If info = 0, b overwritten by the solution vectors, stored columnwise:

if trans = 'N' and mn, rows 1 to n of b contain the least squares solution vectors; the residual sum of squares for the solution in each column is given by the sum of squares of modulus of elements n+1 to m in that column;

if trans = 'N' and m < n, rows 1 to n of b contain the minimum norm solution vectors;

if trans = 'T' or 'C' and mn, rows 1 to m of b contain the minimum norm solution vectors;

if trans = 'T' or 'C' and m < n, rows 1 to m of b contain the least squares solution vectors; the residual sum of squares for the solution in each column is given by the sum of squares of modulus of elements m+1 to n in that column.

Return Values

This function returns a value info.

If info=0, the execution is successful.

If info = -i, the i-th parameter had an illegal value.

If info = i, the i-th diagonal element of the triangular factor of A is zero, so that A does not have full rank; the least squares solution could not be computed.