LAPACK++  2022.07.00
LAPACK C++ API
Standard, AX = B

Solve \(AX \approx B\). More...

Functions

int64_t lapack::gels (lapack::Op trans, int64_t m, int64_t n, int64_t nrhs, double *A, int64_t lda, double *B, int64_t ldb)
 
int64_t lapack::gels (lapack::Op trans, int64_t m, int64_t n, int64_t nrhs, float *A, int64_t lda, float *B, int64_t ldb)
 
int64_t lapack::gels (lapack::Op trans, int64_t m, int64_t n, int64_t nrhs, std::complex< double > *A, int64_t lda, std::complex< double > *B, int64_t ldb)
 Solves overdetermined or underdetermined complex linear systems involving an m-by-n matrix A, or its conjugate-transpose, using a QR or LQ factorization of A. More...
 
int64_t lapack::gels (lapack::Op trans, int64_t m, int64_t n, int64_t nrhs, std::complex< float > *A, int64_t lda, std::complex< float > *B, int64_t ldb)
 
int64_t lapack::gelsd (int64_t m, int64_t n, int64_t nrhs, double *A, int64_t lda, double *B, int64_t ldb, double *S, double rcond, int64_t *rank)
 
int64_t lapack::gelsd (int64_t m, int64_t n, int64_t nrhs, float *A, int64_t lda, float *B, int64_t ldb, float *S, float rcond, int64_t *rank)
 
int64_t lapack::gelsd (int64_t m, int64_t n, int64_t nrhs, std::complex< double > *A, int64_t lda, std::complex< double > *B, int64_t ldb, double *S, double rcond, int64_t *rank)
 Computes the minimum-norm solution to a real linear least squares problem: minimize \(|| b - A x ||_2\) using the singular value decomposition (SVD) of A. More...
 
int64_t lapack::gelsd (int64_t m, int64_t n, int64_t nrhs, std::complex< float > *A, int64_t lda, std::complex< float > *B, int64_t ldb, float *S, float rcond, int64_t *rank)
 
int64_t lapack::gelss (int64_t m, int64_t n, int64_t nrhs, double *A, int64_t lda, double *B, int64_t ldb, double *S, double rcond, int64_t *rank)
 
int64_t lapack::gelss (int64_t m, int64_t n, int64_t nrhs, float *A, int64_t lda, float *B, int64_t ldb, float *S, float rcond, int64_t *rank)
 
int64_t lapack::gelss (int64_t m, int64_t n, int64_t nrhs, std::complex< double > *A, int64_t lda, std::complex< double > *B, int64_t ldb, double *S, double rcond, int64_t *rank)
 Computes the minimum norm solution to a complex linear least squares problem: minimize \(|| b - A x ||_2\) using the singular value decomposition (SVD) of A. More...
 
int64_t lapack::gelss (int64_t m, int64_t n, int64_t nrhs, std::complex< float > *A, int64_t lda, std::complex< float > *B, int64_t ldb, float *S, float rcond, int64_t *rank)
 
int64_t lapack::gelsy (int64_t m, int64_t n, int64_t nrhs, double *A, int64_t lda, double *B, int64_t ldb, int64_t *jpvt, double rcond, int64_t *rank)
 
int64_t lapack::gelsy (int64_t m, int64_t n, int64_t nrhs, float *A, int64_t lda, float *B, int64_t ldb, int64_t *jpvt, float rcond, int64_t *rank)
 
int64_t lapack::gelsy (int64_t m, int64_t n, int64_t nrhs, std::complex< double > *A, int64_t lda, std::complex< double > *B, int64_t ldb, int64_t *jpvt, double rcond, int64_t *rank)
 Computes the minimum-norm solution to a complex linear least squares problem: minimize \(|| A X - B ||_2\) using a complete orthogonal factorization of A. More...
 
int64_t lapack::gelsy (int64_t m, int64_t n, int64_t nrhs, std::complex< float > *A, int64_t lda, std::complex< float > *B, int64_t ldb, int64_t *jpvt, float rcond, int64_t *rank)
 

Detailed Description

Solve \(AX \approx B\).

Function Documentation

◆ gels()

int64_t lapack::gels ( lapack::Op  trans,
int64_t  m,
int64_t  n,
int64_t  nrhs,
std::complex< double > *  A,
int64_t  lda,
std::complex< double > *  B,
int64_t  ldb 
)

Solves overdetermined or underdetermined complex linear systems involving an m-by-n matrix A, or its 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 = NoTrans and m >= n: find the least squares solution of an overdetermined system, i.e., solve the least squares problem minimize \(|| B - A X ||_2\).
  2. If trans = NoTrans and m < n: find the minimum norm solution of an underdetermined system \(A X = B\).
  3. If trans = ConjTrans and m >= n: find the minimum norm solution of an underdetermined system \(A^H X = B\).
  4. If trans = ConjTrans and m < n: find the least squares solution of an overdetermined system, i.e., solve the least squares problem minimize \(|| B - A^H X ||_2\).

Several right hand side vectors b and solution vectors x can be handled in a single call; they are stored as the columns of the m-by-nrhs right hand side matrix B and the n-by-nrhs solution matrix X.

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Parameters
[in]trans
  • lapack::Op::NoTrans: the linear system involves \(A\);
  • lapack::Op::ConjTrans: the linear system involves \(A^H\).
  • lapack::Op::Trans: the linear system involves \(A^T\).
    For real matrices, Trans = ConjTrans. For complex matrices, Trans is illegal.
[in]mThe number of rows of the matrix A. m >= 0.
[in]nThe number of columns of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrices B and X. nrhs >= 0.
[in,out]AThe m-by-n matrix A, stored in an lda-by-n array. On entry, the m-by-n matrix A.
  • If m >= n, A is overwritten by details of its QR factorization as returned by lapack::geqrf;
  • If m < n, A is overwritten by details of its LQ factorization as returned by lapack::gelqf.
[in]ldaThe leading dimension of the array A. lda >= max(1,m).
[in,out]BThe max(m,n)-by-nrhs matrix B, stored in an ldb-by-nrhs array. On entry, the matrix B of right hand side vectors, stored columnwise; B is m-by-nrhs if trans = NoTrans, or n-by-nrhs if trans = ConjTrans. On successful exit, B is overwritten by the solution vectors, stored columnwise:
  • If trans = NoTrans and m >= n, 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 the modulus of elements n+1 to m in that column;
  • If trans = NoTrans and m < n, rows 1 to n of B contain the minimum norm solution vectors;
  • If trans = ConjTrans and m >= n, rows 1 to m of B contain the minimum norm solution vectors;
  • If trans = ConjTrans 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 the modulus of elements m+1 to n in that column.
[in]ldbThe leading dimension of the array B. ldb >= max(1,m,n).
Returns
= 0: successful exit
> 0: if return value = 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.

◆ gelsd()

int64_t lapack::gelsd ( int64_t  m,
int64_t  n,
int64_t  nrhs,
std::complex< double > *  A,
int64_t  lda,
std::complex< double > *  B,
int64_t  ldb,
double *  S,
double  rcond,
int64_t *  rank 
)

Computes the minimum-norm solution to a real linear least squares problem: minimize \(|| b - A x ||_2\) using the singular value decomposition (SVD) of A.

A is an m-by-n matrix which may be rank-deficient.

Several right hand side vectors b and solution vectors x can be handled in a single call; they are stored as the columns of the m-by-nrhs right hand side matrix B and the n-by-nrhs solution matrix X.

The problem is solved in three steps: (1) Reduce the coefficient matrix A to bidiagonal form with Householder transformations, reducing the original problem into a "bidiagonal least squares problem" (BLS) (2) Solve the BLS using a divide and conquer approach. (3) Apply back all the Householder transformations to solve the original least squares problem.

The effective rank of A is determined by treating as zero those singular values which are less than rcond times the largest singular value.

The divide and conquer algorithm makes very mild assumptions about floating point arithmetic. It will work on machines with a guard digit in add/subtract, or on those binary machines without guard digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or Cray-2. It could conceivably fail on hexadecimal or decimal machines without guard digits, but we know of none.

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Parameters
[in]mThe number of rows of the matrix A. m >= 0.
[in]nThe number of columns of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrices B and X. nrhs >= 0.
[in,out]AThe m-by-n matrix A, stored in an lda-by-n array. On entry, the m-by-n matrix A. On exit, A has been destroyed.
[in]ldaThe leading dimension of the array A. lda >= max(1,m).
[in,out]BThe max(m,n)-by-nrhs matrix B or X, stored in an ldb-by-nrhs array. On entry, the m-by-nrhs right hand side matrix B. On exit, B is overwritten by the n-by-nrhs solution matrix X. If m >= n and rank = n, the residual sum-of-squares for the solution in the i-th column is given by the sum of squares of the modulus of elements n+1:m in that column.
[in]ldbThe leading dimension of the array B. ldb >= max(1,m,n).
[out]SThe vector S of length min(m,n). The singular values of A in decreasing order. The condition number of A in the 2-norm = S(1)/S(min(m,n)).
[in]rcondrcond is used to determine the effective rank of A. Singular values S(i) <= rcond*S(1) are treated as zero. If rcond < 0, machine precision is used instead.
[out]rankThe effective rank of A, i.e., the number of singular values which are greater than rcond*S(1).
Returns
= 0: successful exit
> 0: the algorithm for computing the SVD failed to converge; if return value = i, i off-diagonal elements of an intermediate bidiagonal form did not converge to zero.

◆ gelss()

int64_t lapack::gelss ( int64_t  m,
int64_t  n,
int64_t  nrhs,
std::complex< double > *  A,
int64_t  lda,
std::complex< double > *  B,
int64_t  ldb,
double *  S,
double  rcond,
int64_t *  rank 
)

Computes the minimum norm solution to a complex linear least squares problem: minimize \(|| b - A x ||_2\) using the singular value decomposition (SVD) of A.

A is an m-by-n matrix which may be rank-deficient.

Several right hand side vectors b and solution vectors x can be handled in a single call; they are stored as the columns of the m-by-nrhs right hand side matrix B and the n-by-nrhs solution matrix X.

The effective rank of A is determined by treating as zero those singular values which are less than rcond times the largest singular value.

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Parameters
[in]mThe number of rows of the matrix A. m >= 0.
[in]nThe number of columns of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrices B and X. nrhs >= 0.
[in,out]AThe m-by-n matrix A, stored in an lda-by-n array. On entry, the m-by-n matrix A. On exit, the first min(m,n) rows of A are overwritten with its right singular vectors, stored rowwise.
[in]ldaThe leading dimension of the array A. lda >= max(1,m).
[in,out]BThe max(m,n)-by-nrhs matrix B or X, stored in an ldb-by-nrhs array. On entry, the m-by-nrhs right hand side matrix B. On exit, B is overwritten by the n-by-nrhs solution matrix X. If m >= n and rank = n, the residual sum-of-squares for the solution in the i-th column is given by the sum of squares of the modulus of elements n+1:m in that column.
[in]ldbThe leading dimension of the array B. ldb >= max(1,m,n).
[out]SThe vector S of length min(m,n). The singular values of A in decreasing order. The condition number of A in the 2-norm = S(1)/S(min(m,n)).
[in]rcondrcond is used to determine the effective rank of A. Singular values S(i) <= rcond*S(1) are treated as zero. If rcond < 0, machine precision is used instead.
[out]rankThe effective rank of A, i.e., the number of singular values which are greater than rcond*S(1).
Returns
= 0: successful exit
> 0: the algorithm for computing the SVD failed to converge; if return value = i, i off-diagonal elements of an intermediate bidiagonal form did not converge to zero.

◆ gelsy()

int64_t lapack::gelsy ( int64_t  m,
int64_t  n,
int64_t  nrhs,
std::complex< double > *  A,
int64_t  lda,
std::complex< double > *  B,
int64_t  ldb,
int64_t *  jpvt,
double  rcond,
int64_t *  rank 
)

Computes the minimum-norm solution to a complex linear least squares problem: minimize \(|| A X - B ||_2\) using a complete orthogonal factorization of A.

A is an m-by-n matrix which may be rank-deficient.

Several right hand side vectors b and solution vectors x can be handled in a single call; they are stored as the columns of the m-by-nrhs right hand side matrix B and the n-by-nrhs solution matrix X.

The routine first computes a QR factorization with column pivoting:

\[ A P = Q \begin{bmatrix} R_{11} & R_{12} \\ 0 & R_{22} \end{bmatrix}, \]

with R11 defined as the largest leading submatrix whose estimated condition number is less than 1/rcond. The order of R11, rank, is the effective rank of A.

Then, R22 is considered to be negligible, and R12 is annihilated by unitary transformations from the right, arriving at the complete orthogonal factorization:

\[ A P = Q \begin{bmatrix} T_{11} & 0 \\ 0 & 0 \end{bmatrix} Z. \]

The minimum-norm solution is then

\[ X = P Z^H \begin{bmatrix} T_{11}^{-1} Q_1^H B \\ 0 \end{bmatrix}, \]

where \(Q_1\) consists of the first rank columns of Q.

This routine is basically identical to the original gelsx except three differences:

  • The permutation of matrix B (the right hand side) is faster and simpler.
  • The call to the subroutine geqpf has been substituted by the the call to the subroutine geqp3. This subroutine is a BLAS-3 version of the QR factorization with column pivoting.
  • Matrix B (the right hand side) is updated with BLAS-3.

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Parameters
[in]mThe number of rows of the matrix A. m >= 0.
[in]nThe number of columns of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of matrices B and X. nrhs >= 0.
[in,out]AThe m-by-n matrix A, stored in an lda-by-n array. On entry, the m-by-n matrix A. On exit, A has been overwritten by details of its complete orthogonal factorization.
[in]ldaThe leading dimension of the array A. lda >= max(1,m).
[in,out]BThe max(m,n)-by-nrhs matrix B or X, stored in an ldb-by-nrhs array. On entry, the m-by-nrhs right hand side matrix B. On exit, the n-by-nrhs solution matrix X.
[in]ldbThe leading dimension of the array B. ldb >= max(1,m,n).
[in,out]jpvtThe vector jpvt of length n. On entry, if jpvt(i) != 0, the i-th column of A is permuted to the front of AP, otherwise column i is a free column. On exit, if jpvt(i) = k, then the i-th column of \(A P\) was the k-th column of A.
[in]rcondrcond is used to determine the effective rank of A, which is defined as the order of the largest leading triangular submatrix R11 in the QR factorization with pivoting of A, whose estimated condition number < 1/rcond.
[out]rankThe effective rank of A, i.e., the order of the submatrix R11. This is the same as the order of the submatrix T11 in the complete orthogonal factorization of A.
Returns
= 0: successful exit