LAPACK++ 2024.05.31
LAPACK C++ API
|
Functions | |
int64_t | lapack::gtcon (lapack::Norm norm, int64_t n, double const *DL, double const *D, double const *DU, double const *DU2, int64_t const *ipiv, double anorm, double *rcond) |
int64_t | lapack::gtcon (lapack::Norm norm, int64_t n, float const *DL, float const *D, float const *DU, float const *DU2, int64_t const *ipiv, float anorm, float *rcond) |
int64_t | lapack::gtcon (lapack::Norm norm, int64_t n, std::complex< double > const *DL, std::complex< double > const *D, std::complex< double > const *DU, std::complex< double > const *DU2, int64_t const *ipiv, double anorm, double *rcond) |
Estimates the reciprocal of the condition number of a complex tridiagonal matrix A using the LU factorization as computed by lapack::gttrf . | |
int64_t | lapack::gtcon (lapack::Norm norm, int64_t n, std::complex< float > const *DL, std::complex< float > const *D, std::complex< float > const *DU, std::complex< float > const *DU2, int64_t const *ipiv, float anorm, float *rcond) |
int64_t | lapack::gtrfs (lapack::Op trans, int64_t n, int64_t nrhs, double const *DL, double const *D, double const *DU, double const *DLF, double const *DF, double const *DUF, double const *DU2, int64_t const *ipiv, double const *B, int64_t ldb, double *X, int64_t ldx, double *ferr, double *berr) |
int64_t | lapack::gtrfs (lapack::Op trans, int64_t n, int64_t nrhs, float const *DL, float const *D, float const *DU, float const *DLF, float const *DF, float const *DUF, float const *DU2, int64_t const *ipiv, float const *B, int64_t ldb, float *X, int64_t ldx, float *ferr, float *berr) |
int64_t | lapack::gtrfs (lapack::Op trans, int64_t n, int64_t nrhs, std::complex< double > const *DL, std::complex< double > const *D, std::complex< double > const *DU, std::complex< double > const *DLF, std::complex< double > const *DF, std::complex< double > const *DUF, std::complex< double > const *DU2, 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 when the coefficient matrix is tridiagonal, and provides error bounds and backward error estimates for the solution. | |
int64_t | lapack::gtrfs (lapack::Op trans, int64_t n, int64_t nrhs, std::complex< float > const *DL, std::complex< float > const *D, std::complex< float > const *DU, std::complex< float > const *DLF, std::complex< float > const *DF, std::complex< float > const *DUF, std::complex< float > const *DU2, 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::gttrf (int64_t n, double *DL, double *D, double *DU, double *DU2, int64_t *ipiv) |
int64_t | lapack::gttrf (int64_t n, float *DL, float *D, float *DU, float *DU2, int64_t *ipiv) |
int64_t | lapack::gttrf (int64_t n, std::complex< double > *DL, std::complex< double > *D, std::complex< double > *DU, std::complex< double > *DU2, int64_t *ipiv) |
Computes an LU factorization of a complex tridiagonal matrix A using elimination with partial pivoting and row interchanges. | |
int64_t | lapack::gttrf (int64_t n, std::complex< float > *DL, std::complex< float > *D, std::complex< float > *DU, std::complex< float > *DU2, int64_t *ipiv) |
int64_t | lapack::gttrs (lapack::Op trans, int64_t n, int64_t nrhs, double const *DL, double const *D, double const *DU, double const *DU2, int64_t const *ipiv, double *B, int64_t ldb) |
int64_t | lapack::gttrs (lapack::Op trans, int64_t n, int64_t nrhs, float const *DL, float const *D, float const *DU, float const *DU2, int64_t const *ipiv, float *B, int64_t ldb) |
int64_t | lapack::gttrs (lapack::Op trans, int64_t n, int64_t nrhs, std::complex< double > const *DL, std::complex< double > const *D, std::complex< double > const *DU, std::complex< double > const *DU2, int64_t const *ipiv, std::complex< double > *B, int64_t ldb) |
Solves one of the systems of equations. | |
int64_t | lapack::gttrs (lapack::Op trans, int64_t n, int64_t nrhs, std::complex< float > const *DL, std::complex< float > const *D, std::complex< float > const *DU, std::complex< float > const *DU2, int64_t const *ipiv, std::complex< float > *B, int64_t ldb) |
int64_t lapack::gtcon | ( | lapack::Norm | norm, |
int64_t | n, | ||
std::complex< double > const * | DL, | ||
std::complex< double > const * | D, | ||
std::complex< double > const * | DU, | ||
std::complex< double > const * | DU2, | ||
int64_t const * | ipiv, | ||
double | anorm, | ||
double * | rcond | ||
) |
Estimates the reciprocal of the condition number of a complex tridiagonal matrix A using the LU factorization as computed by lapack::gttrf
.
An estimate is obtained for \(|| A^{-1} ||,\) and the reciprocal of the condition number is computed as \(\text{rcond} = 1 / (|| A || \cdot || A^{-1} ||).\)
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[in] | norm | Specifies whether the 1-norm condition number or the infinity-norm condition number is required:
|
[in] | n | The order of the matrix A. n >= 0. |
[in] | DL | The vector DL of length n-1. The (n-1) multipliers that define the matrix L from the LU factorization of A as computed by lapack::gttrf . |
[in] | D | The vector D of length n. The n diagonal elements of the upper triangular matrix U from the LU factorization of A. |
[in] | DU | The vector DU of length n-1. The (n-1) elements of the first superdiagonal of U. |
[in] | DU2 | The vector DU2 of length n-2. The (n-2) elements of the second superdiagonal of U. |
[in] | ipiv | The vector ipiv of length n. The pivot indices; for 1 <= i <= n, row i of the matrix was interchanged with row ipiv(i). ipiv(i) will always be either i or i+1; ipiv(i) = i indicates a row interchange was not required. |
[in] | anorm |
|
[out] | rcond | The reciprocal of the condition number of the matrix A, computed as rcond = 1/(anorm * ainv_norm), where ainv_norm is an estimate of the 1-norm or infinity-norm of \(A^{-1}\) computed in this routine. |
int64_t lapack::gtrfs | ( | lapack::Op | trans, |
int64_t | n, | ||
int64_t | nrhs, | ||
std::complex< double > const * | DL, | ||
std::complex< double > const * | D, | ||
std::complex< double > const * | DU, | ||
std::complex< double > const * | DLF, | ||
std::complex< double > const * | DF, | ||
std::complex< double > const * | DUF, | ||
std::complex< double > const * | DU2, | ||
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 when the coefficient matrix is tridiagonal, 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>
.
[in] | trans | Specifies the form of the system of equations:
|
[in] | n | The order of the matrix A. n >= 0. |
[in] | nrhs | The number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0. |
[in] | DL | The vector DL of length n-1. The (n-1) subdiagonal elements of A. |
[in] | D | The vector D of length n. The diagonal elements of A. |
[in] | DU | The vector DU of length n-1. The (n-1) superdiagonal elements of A. |
[in] | DLF | The vector DLF of length n-1. The (n-1) multipliers that define the matrix L from the LU factorization of A as computed by lapack::gttrf . |
[in] | DF | The vector DF of length n. The n diagonal elements of the upper triangular matrix U from the LU factorization of A. |
[in] | DUF | The vector DUF of length n-1. The (n-1) elements of the first superdiagonal of U. |
[in] | DU2 | The vector DU2 of length n-2. The (n-2) elements of the second superdiagonal of U. |
[in] | ipiv | The vector ipiv of length n. The pivot indices; for 1 <= i <= n, row i of the matrix was interchanged with row ipiv(i). ipiv(i) will always be either i or i+1; ipiv(i) = i indicates a row interchange was not required. |
[in] | B | The n-by-nrhs matrix B, stored in an ldb-by-nrhs array. The right hand side matrix B. |
[in] | ldb | The leading dimension of the array B. ldb >= max(1,n). |
[in,out] | X | The n-by-nrhs matrix X, stored in an ldx-by-nrhs array. On entry, the solution matrix X, as computed by lapack::gttrs . On exit, the improved solution matrix X. |
[in] | ldx | The leading dimension of the array X. ldx >= max(1,n). |
[out] | ferr | The 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] | berr | The 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). |
int64_t lapack::gttrf | ( | int64_t | n, |
std::complex< double > * | DL, | ||
std::complex< double > * | D, | ||
std::complex< double > * | DU, | ||
std::complex< double > * | DU2, | ||
int64_t * | ipiv | ||
) |
Computes an LU factorization of a complex tridiagonal matrix A using elimination with partial pivoting and row interchanges.
The factorization has the form
\[ A = L U \]
where L is a product of permutation and unit lower bidiagonal matrices and U is upper triangular with nonzeros in only the main diagonal and first two superdiagonals.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[in] | n | The order of the matrix A. |
[in,out] | DL | The vector DL of length n-1. On entry, DL must contain the (n-1) sub-diagonal elements of A. On exit, DL is overwritten by the (n-1) multipliers that define the matrix L from the LU factorization of A. |
[in,out] | D | The vector D of length n. On entry, D must contain the diagonal elements of A. On exit, D is overwritten by the n diagonal elements of the upper triangular matrix U from the LU factorization of A. |
[in,out] | DU | The vector DU of length n-1. On entry, DU must contain the (n-1) super-diagonal elements of A. On exit, DU is overwritten by the (n-1) elements of the first super-diagonal of U. |
[out] | DU2 | The vector DU2 of length n-2. On exit, DU2 is overwritten by the (n-2) elements of the second super-diagonal of U. |
[out] | ipiv | The vector ipiv of length n. The pivot indices; for 1 <= i <= n, row i of the matrix was interchanged with row ipiv(i). ipiv(i) will always be either i or i+1; ipiv(i) = i indicates a row interchange was not required. |
int64_t lapack::gttrs | ( | lapack::Op | trans, |
int64_t | n, | ||
int64_t | nrhs, | ||
std::complex< double > const * | DL, | ||
std::complex< double > const * | D, | ||
std::complex< double > const * | DU, | ||
std::complex< double > const * | DU2, | ||
int64_t const * | ipiv, | ||
std::complex< double > * | B, | ||
int64_t | ldb | ||
) |
Solves one of the systems of equations.
\[ A X = B, \]
\[ A^T X = B, \]
or
\[ A^H X = B, \]
with a tridiagonal matrix A using the LU factorization computed by lapack::gttrf
.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[in] | trans | Specifies the form of the system of equations.
|
[in] | n | The order of the matrix A. |
[in] | nrhs | The number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0. |
[in] | DL | The vector DL of length n-1. The (n-1) multipliers that define the matrix L from the LU factorization of A. |
[in] | D | The vector D of length n. The n diagonal elements of the upper triangular matrix U from the LU factorization of A. |
[in] | DU | The vector DU of length n-1. The (n-1) elements of the first super-diagonal of U. |
[in] | DU2 | The vector DU2 of length n-2. The (n-2) elements of the second super-diagonal of U. |
[in] | ipiv | The vector ipiv of length n. The pivot indices; for 1 <= i <= n, row i of the matrix was interchanged with row ipiv(i). ipiv(i) will always be either i or i+1; ipiv(i) = i indicates a row interchange was not required. |
[in,out] | B | The n-by-nrhs matrix B, stored in an ldb-by-nrhs array. On entry, the matrix of right hand side vectors B. On exit, B is overwritten by the solution vectors X. |
[in] | ldb | The leading dimension of the array B. ldb >= max(1,n). |