LAPACK++ 2024.05.31
LAPACK C++ API
Loading...
Searching...
No Matches

Functions

int64_t lapack::geql2 (int64_t m, int64_t n, double *A, int64_t lda, double *tau)
 
int64_t lapack::geql2 (int64_t m, int64_t n, float *A, int64_t lda, float *tau)
 
int64_t lapack::geql2 (int64_t m, int64_t n, std::complex< double > *A, int64_t lda, std::complex< double > *tau)
 Computes a QL factorization of an m-by-n matrix A: \(A = Q L\).
 
int64_t lapack::geql2 (int64_t m, int64_t n, std::complex< float > *A, int64_t lda, std::complex< float > *tau)
 
int64_t lapack::geqlf (int64_t m, int64_t n, double *A, int64_t lda, double *tau)
 
int64_t lapack::geqlf (int64_t m, int64_t n, float *A, int64_t lda, float *tau)
 
int64_t lapack::geqlf (int64_t m, int64_t n, std::complex< double > *A, int64_t lda, std::complex< double > *tau)
 Computes a QL factorization of an m-by-n matrix A: \(A = Q L\).
 
int64_t lapack::geqlf (int64_t m, int64_t n, std::complex< float > *A, int64_t lda, std::complex< float > *tau)
 
int64_t lapack::orgql (int64_t m, int64_t n, int64_t k, double *A, int64_t lda, double const *tau)
 
int64_t lapack::orgql (int64_t m, int64_t n, int64_t k, float *A, int64_t lda, float const *tau)
 
int64_t lapack::ormql (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::ormql (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::ungql (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 last n columns of a product of k elementary reflectors of order m, as returned by lapack::geqlf:
 
int64_t lapack::ungql (int64_t m, int64_t n, int64_t k, std::complex< float > *A, int64_t lda, std::complex< float > const *tau)
 
int64_t lapack::unmql (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::geqlf as follows:
 
int64_t lapack::unmql (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)
 

Detailed Description

Function Documentation

◆ geql2()

int64_t lapack::geql2 ( int64_t  m,
int64_t  n,
std::complex< double > *  A,
int64_t  lda,
std::complex< double > *  tau 
)

Computes a QL factorization of an m-by-n matrix A: \(A = Q L\).

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>.

Parameters
[in]mThe number of rows of the matrix A. m >= 0.
[in]nThe number of columns of the matrix A. n >= 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:
  • if m >= n, the lower triangle of the subarray A(m-n+1:m,1:n) contains the n-by-n lower triangular matrix L;
  • if m <= n, the elements on and below the (n-m)-th superdiagonal contain the m-by-n lower trapezoidal matrix L.
  • The remaining elements, with the array tau, represent the unitary matrix Q as a product of elementary reflectors (see Further Details).
[in]ldaThe leading dimension of the array A. lda >= max(1,m).
[out]tauThe vector tau of length min(m,n). The scalar factors of the elementary reflectors (see Further Details).
Returns
= 0: successful exit
Further Details

The matrix Q is represented as a product of elementary reflectors

\[ Q = H(k) \dots H(2) H(1), \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(m-k+i+1:m) = 0 and v(m-k+i) = 1; v(1:m-k+i-1) is stored on exit in A(1:m-k+i-1,n-k+i), and \(\tau\) in tau(i).

◆ geqlf()

int64_t lapack::geqlf ( int64_t  m,
int64_t  n,
std::complex< double > *  A,
int64_t  lda,
std::complex< double > *  tau 
)

Computes a QL factorization of an m-by-n matrix A: \(A = Q L\).

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>.

Parameters
[in]mThe number of rows of the matrix A. m >= 0.
[in]nThe number of columns of the matrix A. n >= 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:
  • if m >= n, the lower triangle of the subarray A(m-n+1:m,1:n) contains the n-by-n lower triangular matrix L;
  • if m <= n, the elements on and below the (n-m)-th superdiagonal contain the m-by-n lower trapezoidal matrix L.
  • The remaining elements, with the array tau, represent the unitary matrix Q as a product of elementary reflectors (see Further Details).
[in]ldaThe leading dimension of the array A. lda >= max(1,m).
[out]tauThe vector tau of length min(m,n). The scalar factors of the elementary reflectors (see Further Details).
Returns
= 0: successful exit
Further Details

The matrix Q is represented as a product of elementary reflectors

\[ Q = H(k) \dots H(2) H(1), \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(m-k+i+1:m) = 0 and v(m-k+i) = 1; v(1:m-k+i-1) is stored on exit in A(1:m-k+i-1,n-k+i), and \(\tau\) in tau(i).

◆ orgql()

int64_t lapack::orgql ( int64_t  m,
int64_t  n,
int64_t  k,
double *  A,
int64_t  lda,
double const *  tau 
)
See also
lapack::ungql

◆ ormql()

int64_t lapack::ormql ( 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 
)
See also
lapack::unmql

◆ ungql()

int64_t lapack::ungql ( 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 last n columns of a product of k elementary reflectors of order m, as returned by lapack::geqlf:

\[ Q = H(k) \dots H(2) H(1). \]

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>. For real matrices, this is an alias for lapack::orgql.

Parameters
[in]mThe number of rows of the matrix Q. m >= 0.
[in]nThe number of columns of the matrix Q. m >= n >= 0.
[in]kThe number of elementary reflectors whose product defines the matrix Q. n >= k >= 0.
[in,out]AThe m-by-n matrix A, stored in an lda-by-n array. On entry, the (n-k+i)-th column must contain the vector which defines the elementary reflector H(i), for i = 1, 2, ..., k, as returned by lapack::geqlf in the last k columns of its array argument A. On exit, the m-by-n matrix Q.
[in]ldaThe first dimension of the array A. lda >= max(1,m).
[in]tauThe vector tau of length k. tau(i) must contain the scalar factor of the elementary reflector H(i), as returned by lapack::geqlf.
Returns
= 0: successful exit

◆ unmql()

int64_t lapack::unmql ( 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::geqlf as follows:

  • side = Left, trans = NoTrans: \(Q C\)
  • side = Right, trans = NoTrans: \(C Q\)
  • side = Left, trans = ConjTrans: \(Q^H C\)
  • side = Right, trans = ConjTrans: \(C Q^H\)

where Q is a unitary matrix defined as the product of k elementary reflectors, as returned by lapack::geqlf:

\[ Q = H(k) \dots H(2) H(1). \]

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::ormql.

Parameters
[in]side
  • lapack::Side::Left: apply \(Q\) or \(Q^H\) from the Left;
  • lapack::Side::Right: apply \(Q\) or \(Q^H\) from the Right.
[in]trans
  • lapack::Op::NoTrans: No transpose, apply \(Q\);
  • lapack::Op::ConjTrans: Transpose, apply \(Q^H\).
[in]mThe number of rows of the matrix C. m >= 0.
[in]nThe number of columns of the matrix C. n >= 0.
[in]kThe number of elementary reflectors whose product defines the matrix Q.
  • If side = Left, m >= k >= 0;
  • if side = Right, n >= k >= 0.
[in]A
  • If side = Left, the m-by-k matrix A, stored in an lda-by-k array;
  • if side = Right, the n-by-k matrix A, stored in an lda-by-k array.
    The i-th column must contain the vector which defines the elementary reflector H(i), for i = 1, 2, ..., k, as returned by lapack::geqlf in the last k columns of its array argument A.
[in]ldaThe leading dimension of the array A.
  • If side = Left, lda >= max(1,m);
  • if side = Right, lda >= max(1,n).
[in]tauThe vector tau of length k. tau(i) must contain the scalar factor of the elementary reflector H(i), as returned by lapack::geqlf.
[in,out]CThe 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]ldcThe leading dimension of the array C. ldc >= max(1,m).
Returns
= 0: successful exit