LAPACK++ 2024.05.31
LAPACK C++ API
|
Functions | |
int64_t | lapack::geqr (int64_t m, int64_t n, double *A, int64_t lda, double *T, int64_t tsize) |
int64_t | lapack::geqr (int64_t m, int64_t n, float *A, int64_t lda, float *T, int64_t tsize) |
int64_t | lapack::geqr (int64_t m, int64_t n, std::complex< double > *A, int64_t lda, std::complex< double > *T, int64_t tsize) |
Computes a QR factorization of an m-by-n matrix A: | |
int64_t | lapack::geqr (int64_t m, int64_t n, std::complex< float > *A, int64_t lda, std::complex< float > *T, int64_t tsize) |
int64_t | lapack::geqr2 (int64_t m, int64_t n, double *A, int64_t lda, double *tau) |
int64_t | lapack::geqr2 (int64_t m, int64_t n, float *A, int64_t lda, float *tau) |
int64_t | lapack::geqr2 (int64_t m, int64_t n, std::complex< double > *A, int64_t lda, std::complex< double > *tau) |
Computes a QR factorization of an m-by-n matrix A: | |
int64_t | lapack::geqr2 (int64_t m, int64_t n, std::complex< float > *A, int64_t lda, std::complex< float > *tau) |
int64_t | lapack::geqrf (int64_t m, int64_t n, double *A, int64_t lda, double *tau) |
int64_t | lapack::geqrf (int64_t m, int64_t n, float *A, int64_t lda, float *tau) |
int64_t | lapack::geqrf (int64_t m, int64_t n, std::complex< double > *A, int64_t lda, std::complex< double > *tau) |
Computes a QR factorization of an m-by-n matrix A: \(A = Q R\). | |
int64_t | lapack::geqrf (int64_t m, int64_t n, std::complex< float > *A, int64_t lda, std::complex< float > *tau) |
int64_t | lapack::orgqr (int64_t m, int64_t n, int64_t k, double *A, int64_t lda, double const *tau) |
int64_t | lapack::orgqr (int64_t m, int64_t n, int64_t k, float *A, int64_t lda, float const *tau) |
int64_t | lapack::ormqr (lapack::Side side, lapack::Op trans, int64_t m, int64_t n, int64_t k, double const *A, int64_t lda, double const *tau, double *C, int64_t ldc) |
int64_t | lapack::ormqr (lapack::Side side, lapack::Op trans, int64_t m, int64_t n, int64_t k, float const *A, int64_t lda, float const *tau, float *C, int64_t ldc) |
int64_t | lapack::ungqr (int64_t m, int64_t n, int64_t k, std::complex< double > *A, int64_t lda, std::complex< double > const *tau) |
Generates an m-by-n matrix Q with orthonormal columns, which is defined as the first n columns of a product of k elementary reflectors of order m, as returned by lapack::geqrf : | |
int64_t | lapack::ungqr (int64_t m, int64_t n, int64_t k, std::complex< float > *A, int64_t lda, std::complex< float > const *tau) |
int64_t | lapack::unmqr (lapack::Side side, lapack::Op trans, int64_t m, int64_t n, int64_t k, std::complex< double > const *A, int64_t lda, std::complex< double > const *tau, std::complex< double > *C, int64_t ldc) |
Multiplies the general m-by-n matrix C by Q from lapack::geqrf as follows: | |
int64_t | lapack::unmqr (lapack::Side side, lapack::Op trans, int64_t m, int64_t n, int64_t k, std::complex< float > const *A, int64_t lda, std::complex< float > const *tau, std::complex< float > *C, int64_t ldc) |
int64_t lapack::geqr | ( | int64_t | m, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
std::complex< double > * | T, | ||
int64_t | tsize | ||
) |
Computes a QR factorization of an m-by-n matrix A:
\[ A = Q \begin{bmatrix} R \\ 0 \end{bmatrix}, \]
Q is a m-by-m orthogonal matrix; R is an upper-triangular n-by-n matrix; 0 is a (m - n)-by-n zero matrix, if m > n.
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,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 elements on and above the diagonal of the array contain the min(m,n)-by-n upper trapezoidal matrix R (R is upper triangular if m >= n); the elements below the diagonal are used to store part of the data structure to represent Q. |
[in] | lda | The leading dimension of the array A. lda >= max(1,m). |
[out] | T | The vector T of length max(5,tsize). On successful exit, T[0] returns optimal (or either minimal or optimal, if query is assumed) tsize. See tsize for details. Remaining T contains part of the data structure used to represent Q. If one wants to apply or construct Q, then one needs to keep T (in addition to A) and pass it to further subroutines. |
[in] | tsize | If tsize >= 5, the dimension of the array T. If tsize = -1 or -2, then a workspace query is assumed. The routine only calculates the sizes of the T array, returns this value as the first entries of the T array, and no error message related to T is issued. If tsize = -1, the routine calculates optimal size of T for the optimum performance and returns this value in T[0]. If tsize = -2, the routine calculates minimal size of T and returns this value in T[0]. |
= | 0: successful exit |
The goal of the interface is to give maximum freedom to the developers for creating any QR factorization algorithm they wish. The trapezoidal R has to be stored in the upper part of A. The lower part of A and the array T can be used to store any relevant information for applying or constructing the Q factor.
Caution: One should not expect the size of T to be the same from one LAPACK implementation to the other, or even from one execution to the other. A workspace query for T is needed at each execution. However, for a given execution, the size of T are fixed and will not change from one query to the next.
These details are particular for the Netlib LAPACK implementation. Users should not take them for granted. These details may change in the future, and are not likely true for another LAPACK implementation. These details are relevant if one wants to try to understand the code. They are not part of the interface.
In this version,
T[1]: row block size (mb) T[2]: column block size (nb) T[5:TSIZE-1]: data structure needed for Q, computed by latsqr or geqrt
Depending on the matrix dimensions m and n, and row and column block sizes mb and nb returned by ilaenv, geqr will use either latsqr (if the matrix is tall-and-skinny) or geqrt to compute the QR factorization.
int64_t lapack::geqr2 | ( | int64_t | m, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
std::complex< double > * | tau | ||
) |
Computes a QR factorization of an m-by-n matrix A:
\[ A = Q R. \]
This is the unblocked Level 2 BLAS version of the algorithm.
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,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 elements on and above the diagonal of the array contain the min(m,n)-by-n upper trapezoidal matrix R (R is upper triangular if m >= n). The elements below the diagonal, with the array tau, represent the unitary matrix Q as a product of min(m,n) elementary reflectors (see Further Details). |
[in] | lda | The leading dimension of the array A. lda >= max(1,m). |
[out] | tau | The vector tau of length min(m,n). The scalar factors of the elementary reflectors (see Further Details). |
The matrix Q is represented as a product of elementary reflectors
\[ Q = H(1) H(2) \dots H(k) \text{ where } k = \min(m,n). \]
Each H(i) has the form
\[ H(i) = I - \tau v v^H \]
where \(\tau\) is a scalar, and v is a vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and \(\tau\) in tau(i).
int64_t lapack::geqrf | ( | int64_t | m, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
std::complex< double > * | tau | ||
) |
Computes a QR factorization of an m-by-n matrix A: \(A = Q R\).
This is the blocked Level 3 BLAS version of the algorithm.
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,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 elements on and above the diagonal of the array contain the min(m,n)-by-n upper trapezoidal matrix R (R is upper triangular if m >= n). The elements below the diagonal, with the array tau, represent the unitary matrix Q as a product of min(m,n) elementary reflectors (see Further Details). |
[in] | lda | The leading dimension of the array A. lda >= max(1,m). |
[out] | tau | The vector tau of length min(m,n). The scalar factors of the elementary reflectors (see Further Details). |
The matrix Q is represented as a product of elementary reflectors
\[ Q = H(1) H(2) \dots H(k) \text{ where } k = \min(m,n). \]
Each H(i) has the form
\[ H(i) = I - \tau v v^H \]
where \(\tau\) is a scalar, and v is a vector with v(1:i-1) = 0 and v(i) = 1; v(i+1:m) is stored on exit in A(i+1:m,i), and \(\tau\) in tau(i).
int64_t lapack::orgqr | ( | int64_t | m, |
int64_t | n, | ||
int64_t | k, | ||
double * | A, | ||
int64_t | lda, | ||
double const * | tau | ||
) |
int64_t lapack::ormqr | ( | lapack::Side | side, |
lapack::Op | trans, | ||
int64_t | m, | ||
int64_t | n, | ||
int64_t | k, | ||
double const * | A, | ||
int64_t | lda, | ||
double const * | tau, | ||
double * | C, | ||
int64_t | ldc | ||
) |
int64_t lapack::ungqr | ( | int64_t | m, |
int64_t | n, | ||
int64_t | k, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
std::complex< double > const * | tau | ||
) |
Generates an m-by-n matrix Q with orthonormal columns, which is defined as the first n columns of a product of k elementary reflectors of order m, as returned by lapack::geqrf
:
\[ Q = H(1) H(2) \dots H(k). \]
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this is an alias for lapack::orgqr
.
[in] | m | The number of rows of the matrix Q. m >= 0. |
[in] | n | The number of columns of the matrix Q. m >= n >= 0. |
[in] | k | The number of elementary reflectors whose product defines the matrix Q. n >= k >= 0. |
[in,out] | A | The m-by-n matrix A, stored in an lda-by-n array. On entry, the i-th column must contain the vector which defines the elementary reflector H(i), for i = 1, 2, ..., k, as returned by lapack::geqrf in the first k columns of its array argument A. On exit, the m-by-n matrix Q. |
[in] | lda | The first dimension of the array A. lda >= max(1,m). |
[in] | tau | The vector tau of length k. tau(i) must contain the scalar factor of the elementary reflector H(i), as returned by lapack::geqrf . |
int64_t lapack::unmqr | ( | lapack::Side | side, |
lapack::Op | trans, | ||
int64_t | m, | ||
int64_t | n, | ||
int64_t | k, | ||
std::complex< double > const * | A, | ||
int64_t | lda, | ||
std::complex< double > const * | tau, | ||
std::complex< double > * | C, | ||
int64_t | ldc | ||
) |
Multiplies the general m-by-n matrix C by Q from lapack::geqrf
as follows:
where Q is a unitary matrix defined as the product of k elementary reflectors, as returned by lapack::geqrf
:
\[ Q = H(1) H(2) \dots H(k). \]
Q is of order m if side = Left and of order n if side = Right.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this is an alias for lapack::ormqr
.
[in] | side |
|
[in] | trans |
|
[in] | m | The number of rows of the matrix C. m >= 0. |
[in] | n | The number of columns of the matrix C. n >= 0. |
[in] | k | The number of elementary reflectors whose product defines the matrix Q.
|
[in] | A |
|
[in] | lda | The leading dimension of the array A.
|
[in] | tau | The vector tau of length k. tau(i) must contain the scalar factor of the elementary reflector H(i), as returned by lapack::geqrf . |
[in,out] | C | The m-by-n matrix C, stored in an ldc-by-n array. On entry, the m-by-n matrix C. On exit, C is overwritten by \(Q C\) or \(Q^H C\) or \(C Q^H\) or \(C Q\). |
[in] | ldc | The leading dimension of the array C. ldc >= max(1,m). |