LAPACK++
2022.05.00
LAPACK C++ API
|
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) |
Solve \(AX \approx B\).
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:
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>
.
[in] | trans |
|
[in] | m | The number of rows of the matrix A. m >= 0. |
[in] | n | The number of columns of the matrix A. n >= 0. |
[in] | nrhs | The number of right hand sides, i.e., the number of columns of the matrices B and X. nrhs >= 0. |
[in,out] | A | The m-by-n matrix A, stored in an lda-by-n array. On entry, the m-by-n matrix A.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,m). |
[in,out] | B | The 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:
|
[in] | ldb | The leading dimension of the array B. ldb >= max(1,m,n). |
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>
.
[in] | m | The number of rows of the matrix A. m >= 0. |
[in] | n | The number of columns of the matrix A. n >= 0. |
[in] | nrhs | The number of right hand sides, i.e., the number of columns of the matrices B and X. nrhs >= 0. |
[in,out] | A | The 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] | lda | The leading dimension of the array A. lda >= max(1,m). |
[in,out] | B | The 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] | ldb | The leading dimension of the array B. ldb >= max(1,m,n). |
[out] | S | The 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] | rcond | rcond 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] | rank | The effective rank of A, i.e., the number of singular values which are greater than rcond*S(1). |
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>
.
[in] | m | The number of rows of the matrix A. m >= 0. |
[in] | n | The number of columns of the matrix A. n >= 0. |
[in] | nrhs | The number of right hand sides, i.e., the number of columns of the matrices B and X. nrhs >= 0. |
[in,out] | A | The 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] | lda | The leading dimension of the array A. lda >= max(1,m). |
[in,out] | B | The 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] | ldb | The leading dimension of the array B. ldb >= max(1,m,n). |
[out] | S | The 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] | rcond | rcond 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] | rank | The effective rank of A, i.e., the number of singular values which are greater than rcond*S(1). |
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:
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[in] | m | The number of rows of the matrix A. m >= 0. |
[in] | n | The number of columns of the matrix A. n >= 0. |
[in] | nrhs | The number of right hand sides, i.e., the number of columns of matrices B and X. nrhs >= 0. |
[in,out] | A | The 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] | lda | The leading dimension of the array A. lda >= max(1,m). |
[in,out] | B | The 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] | ldb | The leading dimension of the array B. ldb >= max(1,m,n). |
[in,out] | jpvt | The 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] | rcond | rcond 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] | rank | The 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. |