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

Functions

int64_t lapack::gecon (lapack::Norm norm, int64_t n, double const *A, int64_t lda, double anorm, double *rcond)
 
int64_t lapack::gecon (lapack::Norm norm, int64_t n, float const *A, int64_t lda, float anorm, float *rcond)
 
int64_t lapack::gecon (lapack::Norm norm, int64_t n, std::complex< double > const *A, int64_t lda, double anorm, double *rcond)
 Estimates the reciprocal of the condition number of a general matrix A, in either the 1-norm or the infinity-norm, using the LU factorization computed by lapack::getrf.
 
int64_t lapack::gecon (lapack::Norm norm, int64_t n, std::complex< float > const *A, int64_t lda, float anorm, float *rcond)
 
int64_t lapack::geequ (int64_t m, int64_t n, double const *A, int64_t lda, double *R, double *C, double *rowcnd, double *colcnd, double *amax)
 
int64_t lapack::geequ (int64_t m, int64_t n, float const *A, int64_t lda, float *R, float *C, float *rowcnd, float *colcnd, float *amax)
 
int64_t lapack::geequ (int64_t m, int64_t n, std::complex< double > const *A, int64_t lda, double *R, double *C, double *rowcnd, double *colcnd, double *amax)
 Computes row and column scalings intended to equilibrate an m-by-n matrix A and reduce its condition number.
 
int64_t lapack::geequ (int64_t m, int64_t n, std::complex< float > const *A, int64_t lda, float *R, float *C, float *rowcnd, float *colcnd, float *amax)
 
int64_t lapack::geequb (int64_t m, int64_t n, double const *A, int64_t lda, double *R, double *C, double *rowcnd, double *colcnd, double *amax)
 
int64_t lapack::geequb (int64_t m, int64_t n, float const *A, int64_t lda, float *R, float *C, float *rowcnd, float *colcnd, float *amax)
 
int64_t lapack::geequb (int64_t m, int64_t n, std::complex< double > const *A, int64_t lda, double *R, double *C, double *rowcnd, double *colcnd, double *amax)
 Computes row and column scalings intended to equilibrate an m-by-n matrix A and reduce its condition number.
 
int64_t lapack::geequb (int64_t m, int64_t n, std::complex< float > const *A, int64_t lda, float *R, float *C, float *rowcnd, float *colcnd, float *amax)
 
int64_t lapack::gerfs (lapack::Op trans, int64_t n, int64_t nrhs, double const *A, int64_t lda, double const *AF, int64_t ldaf, int64_t const *ipiv, double const *B, int64_t ldb, double *X, int64_t ldx, double *ferr, double *berr)
 
int64_t lapack::gerfs (lapack::Op trans, int64_t n, int64_t nrhs, float const *A, int64_t lda, float const *AF, int64_t ldaf, int64_t const *ipiv, float const *B, int64_t ldb, float *X, int64_t ldx, float *ferr, float *berr)
 
int64_t lapack::gerfs (lapack::Op trans, int64_t n, int64_t nrhs, std::complex< double > const *A, int64_t lda, std::complex< double > const *AF, int64_t ldaf, int64_t const *ipiv, std::complex< double > const *B, int64_t ldb, std::complex< double > *X, int64_t ldx, double *ferr, double *berr)
 Improves the computed solution to a system of linear equations and provides error bounds and backward error estimates for the solution.
 
int64_t lapack::gerfs (lapack::Op trans, int64_t n, int64_t nrhs, std::complex< float > const *A, int64_t lda, std::complex< float > const *AF, int64_t ldaf, int64_t const *ipiv, std::complex< float > const *B, int64_t ldb, std::complex< float > *X, int64_t ldx, float *ferr, float *berr)
 
int64_t lapack::getf2 (int64_t m, int64_t n, double *A, int64_t lda, int64_t *ipiv)
 
int64_t lapack::getf2 (int64_t m, int64_t n, float *A, int64_t lda, int64_t *ipiv)
 
int64_t lapack::getf2 (int64_t m, int64_t n, std::complex< double > *A, int64_t lda, int64_t *ipiv)
 Computes an LU factorization of a general m-by-n matrix A using partial pivoting with row interchanges.
 
int64_t lapack::getf2 (int64_t m, int64_t n, std::complex< float > *A, int64_t lda, int64_t *ipiv)
 
int64_t lapack::getrf (int64_t m, int64_t n, double *A, int64_t lda, int64_t *ipiv)
 
int64_t lapack::getrf (int64_t m, int64_t n, float *A, int64_t lda, int64_t *ipiv)
 
int64_t lapack::getrf (int64_t m, int64_t n, std::complex< double > *A, int64_t lda, int64_t *ipiv)
 Computes an LU factorization of a general m-by-n matrix A using partial pivoting with row interchanges.
 
int64_t lapack::getrf (int64_t m, int64_t n, std::complex< float > *A, int64_t lda, int64_t *ipiv)
 
int64_t lapack::getrf2 (int64_t m, int64_t n, double *A, int64_t lda, int64_t *ipiv)
 
int64_t lapack::getrf2 (int64_t m, int64_t n, float *A, int64_t lda, int64_t *ipiv)
 
int64_t lapack::getrf2 (int64_t m, int64_t n, std::complex< double > *A, int64_t lda, int64_t *ipiv)
 Computes an LU factorization of a general m-by-n matrix A using partial pivoting with row interchanges.
 
int64_t lapack::getrf2 (int64_t m, int64_t n, std::complex< float > *A, int64_t lda, int64_t *ipiv)
 
int64_t lapack::getri (int64_t n, double *A, int64_t lda, int64_t const *ipiv)
 
int64_t lapack::getri (int64_t n, float *A, int64_t lda, int64_t const *ipiv)
 
int64_t lapack::getri (int64_t n, std::complex< double > *A, int64_t lda, int64_t const *ipiv)
 Computes the inverse of a matrix using the LU factorization computed by lapack::getrf.
 
int64_t lapack::getri (int64_t n, std::complex< float > *A, int64_t lda, int64_t const *ipiv)
 
int64_t lapack::getrs (lapack::Op trans, int64_t n, int64_t nrhs, double const *A, int64_t lda, int64_t const *ipiv, double *B, int64_t ldb)
 
int64_t lapack::getrs (lapack::Op trans, int64_t n, int64_t nrhs, float const *A, int64_t lda, int64_t const *ipiv, float *B, int64_t ldb)
 
int64_t lapack::getrs (lapack::Op trans, int64_t n, int64_t nrhs, std::complex< double > const *A, int64_t lda, int64_t const *ipiv, std::complex< double > *B, int64_t ldb)
 Solves a system of linear equations.
 
int64_t lapack::getrs (lapack::Op trans, int64_t n, int64_t nrhs, std::complex< float > const *A, int64_t lda, int64_t const *ipiv, std::complex< float > *B, int64_t ldb)
 
void lapack::laswp (int64_t n, double *A, int64_t lda, int64_t k1, int64_t k2, int64_t const *ipiv, int64_t incx)
 
void lapack::laswp (int64_t n, float *A, int64_t lda, int64_t k1, int64_t k2, int64_t const *ipiv, int64_t incx)
 
void lapack::laswp (int64_t n, std::complex< double > *A, int64_t lda, int64_t k1, int64_t k2, int64_t const *ipiv, int64_t incx)
 Performs a series of row interchanges on the matrix A.
 
void lapack::laswp (int64_t n, std::complex< float > *A, int64_t lda, int64_t k1, int64_t k2, int64_t const *ipiv, int64_t incx)
 

Detailed Description

Function Documentation

◆ gecon()

int64_t lapack::gecon ( lapack::Norm  norm,
int64_t  n,
std::complex< double > const *  A,
int64_t  lda,
double  anorm,
double *  rcond 
)

Estimates the reciprocal of the condition number of a general matrix A, in either the 1-norm or the infinity-norm, using the LU factorization computed by lapack::getrf.

An estimate is obtained for norm(inv(A)), and the reciprocal of the condition number is computed as rcond = 1 / ( norm(A) * norm(inv(A)) ).

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

Parameters
[in]normWhether the 1-norm condition number or the infinity-norm condition number is required:
  • lapack::Norm::One: 1-norm;
  • lapack::Norm::Inf: Infinity-norm.
[in]nThe order of the matrix A. n >= 0.
[in]AThe n-by-n matrix A, stored in an lda-by-n array. The factors L and U from the factorization \(A = P L U\) as computed by lapack::getrf.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in]anorm
  • If norm = One, the 1-norm of the original matrix A.
  • If norm = Inf, the infinity-norm of the original matrix A.
[out]rcondThe reciprocal of the condition number of the matrix A, computed as rcond = 1/(norm(A) * norm(inv(A))).
Returns
= 0: successful exit

◆ geequ()

int64_t lapack::geequ ( int64_t  m,
int64_t  n,
std::complex< double > const *  A,
int64_t  lda,
double *  R,
double *  C,
double *  rowcnd,
double *  colcnd,
double *  amax 
)

Computes row and column scalings intended to equilibrate an m-by-n matrix A and reduce its condition number.

R returns the row scale factors and C the column scale factors, chosen to try to make the largest element in each row and column of the matrix B with elements \(B_{i,j} = R_{i} A_{i,j} C_{j}\) have absolute value 1.

\(R_{i}\) and \(C_{j}\) are restricted to be between smlnum = smallest safe number and bignum = largest safe number. Use of these scaling factors is not guaranteed to reduce the condition number of A but works well in practice.

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]AThe m-by-n matrix A, stored in an lda-by-n array. The m-by-n matrix whose equilibration factors are to be computed.
[in]ldaThe leading dimension of the array A. lda >= max(1,m).
[out]RThe vector R of length m. If successful or return value > m, R contains the row scale factors for A.
[out]CThe vector C of length n. If successful, C contains the column scale factors for A.
[out]rowcndIf successful or return value > m, rowcnd contains the ratio of the smallest R(i) to the largest R(i). If rowcnd >= 0.1 and amax is neither too large nor too small, it is not worth scaling by R.
[out]colcndIf successful, colcnd contains the ratio of the smallest C(i) to the largest C(i). If colcnd >= 0.1, it is not worth scaling by C.
[out]amaxAbsolute value of largest matrix element. If amax is very close to overflow or very close to underflow, the matrix should be scaled.
Returns
= 0: successful exit
> 0 and <= m: if return value = i, the i-th row of A is exactly zero
> m: if return value = i, the (i-m)-th column of A is exactly zero

◆ geequb()

int64_t lapack::geequb ( int64_t  m,
int64_t  n,
std::complex< double > const *  A,
int64_t  lda,
double *  R,
double *  C,
double *  rowcnd,
double *  colcnd,
double *  amax 
)

Computes row and column scalings intended to equilibrate an m-by-n matrix A and reduce its condition number.

R returns the row scale factors and C the column scale factors, chosen to try to make the largest element in each row and column of the matrix B with elements \(B_{i,j} = R_{i} A_{i,j} C_{j}\) have an absolute value of at most the radix.

\(R_{i}\) and \(C_{j}\) are restricted to be a power of the radix between smlnum = smallest safe number and bignum = largest safe number. Use of these scaling factors is not guaranteed to reduce the condition number of A but works well in practice.

This routine differs from lapack::geequ by restricting the scaling factors to a power of the radix. Barring over- and underflow, scaling by these factors introduces no additional rounding errors. However, the scaled entries' magnitudes are no longer approximately 1 but lie between sqrt(radix) and 1/sqrt(radix).

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]AThe m-by-n matrix A, stored in an lda-by-n array. The m-by-n matrix whose equilibration factors are to be computed.
[in]ldaThe leading dimension of the array A. lda >= max(1,m).
[out]RThe vector R of length m. If successful or return value > m, R contains the row scale factors for A.
[out]CThe vector C of length n. If successful, C contains the column scale factors for A.
[out]rowcndIf successful or return value > m, rowcnd contains the ratio of the smallest R(i) to the largest R(i). If rowcnd >= 0.1 and amax is neither too large nor too small, it is not worth scaling by R.
[out]colcndIf successful, colcnd contains the ratio of the smallest C(i) to the largest C(i). If colcnd >= 0.1, it is not worth scaling by C.
[out]amaxAbsolute value of largest matrix element. If amax is very close to overflow or very close to underflow, the matrix should be scaled.
Returns
= 0: successful exit
> 0 and <= m: if return value = i, the i-th row of A is exactly zero
> m: if return value = i, the (i-m)-th column of A is exactly zero

◆ gerfs()

int64_t lapack::gerfs ( lapack::Op  trans,
int64_t  n,
int64_t  nrhs,
std::complex< double > const *  A,
int64_t  lda,
std::complex< double > const *  AF,
int64_t  ldaf,
int64_t const *  ipiv,
std::complex< double > const *  B,
int64_t  ldb,
std::complex< double > *  X,
int64_t  ldx,
double *  ferr,
double *  berr 
)

Improves the computed solution to a system of linear equations and provides error bounds and backward error estimates for the solution.

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

Parameters
[in]transThe form of the system of equations:
  • lapack::Op::NoTrans: \(A X = B\) (No transpose)
  • lapack::Op::Trans: \(A^T X = B\) (Transpose)
  • lapack::Op::ConjTrans: \(A^H X = B\) (Conjugate transpose)
[in]nThe order 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]AThe n-by-n matrix A, stored in an lda-by-n array. The original n-by-n matrix A.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in]AFThe n-by-n matrix AF, stored in an ldaf-by-n array. The factors L and U from the factorization \(A = P L U\) as computed by lapack::getrf.
[in]ldafThe leading dimension of the array AF. ldaf >= max(1,n).
[in]ipivThe vector ipiv of length n. The pivot indices from lapack::getrf; for 1 <= i <= n, row i of the matrix was interchanged with row ipiv(i).
[in]BThe n-by-nrhs matrix B, stored in an ldb-by-nrhs array. The right hand side matrix B.
[in]ldbThe leading dimension of the array B. ldb >= max(1,n).
[in,out]XThe n-by-nrhs matrix X, stored in an ldx-by-nrhs array. On entry, the solution matrix X, as computed by lapack::getrs. On exit, the improved solution matrix X.
[in]ldxThe leading dimension of the array X. ldx >= max(1,n).
[out]ferrThe vector ferr of length nrhs. The estimated forward error bound for each solution vector X(j) (the j-th column of the solution matrix X). If XTRUE is the true solution corresponding to X(j), ferr(j) is an estimated upper bound for the magnitude of the largest element in (X(j) - XTRUE) divided by the magnitude of the largest element in X(j). The estimate is as reliable as the estimate for rcond, and is almost always a slight overestimate of the true error.
[out]berrThe vector berr of length nrhs. The componentwise relative backward error of each solution vector X(j) (i.e., the smallest relative change in any element of A or B that makes X(j) an exact solution).
Returns
= 0: successful exit

◆ getf2()

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

Computes an LU factorization of a general m-by-n matrix A using partial pivoting with row interchanges.

The factorization has the form

\[ A = P L U \]

where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the right-looking 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 to be factored. On exit, the factors L and U from the factorization \(A = P L U;\) the unit diagonal elements of L are not stored.
[in]ldaThe leading dimension of the array A. lda >= max(1,m).
[out]ipivThe vector ipiv of length min(m,n). The pivot indices; for 1 <= i <= min(m,n), row i of the matrix was interchanged with row ipiv(i).
Returns
= 0: successful exit
> 0: if return value = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.

◆ getrf()

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

Computes an LU factorization of a general m-by-n matrix A using partial pivoting with row interchanges.

The factorization has the form

\[ A = P L U \]

where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the right-looking 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 to be factored. On exit, the factors L and U from the factorization \(A = P L U;\) the unit diagonal elements of L are not stored.
[in]ldaThe leading dimension of the array A. lda >= max(1,m).
[out]ipivThe vector ipiv of length min(m,n). The pivot indices; for 1 <= i <= min(m,n), row i of the matrix was interchanged with row ipiv(i).
Returns
= 0: successful exit
> 0: if return value = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.

◆ getrf2()

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

Computes an LU factorization of a general m-by-n matrix A using partial pivoting with row interchanges.

The factorization has the form

\[ A = P L U \]

where P is a permutation matrix, L is lower triangular with unit diagonal elements (lower trapezoidal if m > n), and U is upper triangular (upper trapezoidal if m < n).

This is the recursive version of the algorithm. It divides the matrix into four submatrices:

\[ A = \begin{bmatrix} A_{11} & A_{12} \\ A_{21} & A_{22} \end{bmatrix} \]

where \(A_{11}\) is n1-by-n1 and \(A_{22}\) is n2-by-n2, with n1 = min(m,n)/2 and n2 = n-n1. The subroutine calls itself to factor

\[ \begin{bmatrix} A_{11} \\ A_{21} \end{bmatrix}, \]

does the swaps on

\[ \begin{bmatrix} A_{12} \\ A_{22} \end{bmatrix}, \]

solves \(A_{12},\) updates \(A_{22},\) calls itself to factor \(A_{22},\) and does the swaps on \(A_{21}.\)

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 to be factored. On exit, the factors L and U from the factorization \(A = P L U;\) the unit diagonal elements of L are not stored.
[in]ldaThe leading dimension of the array A. lda >= max(1,m).
[out]ipivThe vector ipiv of length min(m,n). The pivot indices; for 1 <= i <= min(m,n), row i of the matrix was interchanged with row ipiv(i).
Returns
= 0: successful exit
> 0: if return value = i, U(i,i) is exactly zero. The factorization has been completed, but the factor U is exactly singular, and division by zero will occur if it is used to solve a system of equations.

◆ getri()

int64_t lapack::getri ( int64_t  n,
std::complex< double > *  A,
int64_t  lda,
int64_t const *  ipiv 
)

Computes the inverse of a matrix using the LU factorization computed by lapack::getrf.

This method inverts U and then computes \(A^{-1}\) by solving the system

\[ A^{-1} L = U^{-1} \text{ for } A^{-1}. \]

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

Parameters
[in]nThe order of the matrix A. n >= 0.
[in,out]AThe n-by-n matrix A, stored in an lda-by-n array. On entry, the factors L and U from the factorization \(A = P L U\) as computed by lapack::getrf. On successful exit, the inverse of the original matrix A.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in]ipivThe vector ipiv of length n. The pivot indices from lapack::getrf; for 1 <= i <= n, row i of the matrix was interchanged with row ipiv(i).
Returns
= 0: successful exit
> 0: if return value = i, U(i,i) is exactly zero; the matrix is singular and its inverse could not be computed.

◆ getrs()

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

Solves a system of linear equations.

\[ A X = B, \]

\[ A^T X = B, \]

or

\[ A^H X = B \]

with a general n-by-n matrix A using the LU factorization computed by lapack::getrf.

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

Parameters
[in]transThe form of the system of equations:
  • lapack::Op::NoTrans: \(A X = B\) (No transpose)
  • lapack::Op::Trans: \(A^T X = B\) (Transpose)
  • lapack::Op::ConjTrans: \(A^H X = B\) (Conjugate transpose)
[in]nThe order of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0.
[in]AThe n-by-n matrix A, stored in an lda-by-n array. The factors L and U from the factorization \(A = P L U\) as computed by lapack::getrf.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in]ipivThe vector ipiv of length n. The pivot indices from lapack::getrf; for 1 <= i <= n, row i of the matrix was interchanged with row ipiv(i).
[in,out]BThe n-by-nrhs matrix B, stored in an ldb-by-nrhs array. On entry, the right hand side matrix B. On exit, the solution matrix X.
[in]ldbThe leading dimension of the array B. ldb >= max(1,n).
Returns
= 0: successful exit

◆ laswp()

void lapack::laswp ( int64_t  n,
std::complex< double > *  A,
int64_t  lda,
int64_t  k1,
int64_t  k2,
int64_t const *  ipiv,
int64_t  incx 
)

Performs a series of row interchanges on the matrix A.

One row interchange is initiated for each of rows k1 through k2 of A.

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

Parameters
[in]nThe number of columns of the matrix A.
[in,out]AThe m-by-n matrix A, stored in an lda-by-n array. On entry, the matrix of column dimension n to which the row interchanges will be applied. On exit, the permuted matrix. Note that the number of rows, m, is implicit in ipiv; m <= lda.
[in]ldaThe leading dimension of the array A.
[in]k1The first element of ipiv for which a row interchange will be done.
[in]k2(k2-k1+1) is the number of elements of ipiv for which a row interchange will be done.
[in]ipivThe vector ipiv of length k1+(k2-k1)*abs(incx). The vector of pivot indices. Only the elements in positions k1 through k1+(k2-k1)*abs(incx) of ipiv are accessed. ipiv(k1+( \(K-\)k1)*abs(incx)) = L implies rows K and L are to be interchanged.
[in]incxThe increment between successive values of ipiv. If incx is negative, the pivots are applied in reverse order.