LAPACK++ 2024.05.31
LAPACK C++ API
|
Functions | |
int64_t | lapack::gbcon (lapack::Norm norm, int64_t n, int64_t kl, int64_t ku, double const *AB, int64_t ldab, int64_t const *ipiv, double anorm, double *rcond) |
int64_t | lapack::gbcon (lapack::Norm norm, int64_t n, int64_t kl, int64_t ku, float const *AB, int64_t ldab, int64_t const *ipiv, float anorm, float *rcond) |
int64_t | lapack::gbcon (lapack::Norm norm, int64_t n, int64_t kl, int64_t ku, std::complex< double > const *AB, int64_t ldab, int64_t const *ipiv, double anorm, double *rcond) |
Estimates the reciprocal of the condition number of a general band matrix A, in either the 1-norm or the infinity-norm, using the LU factorization computed by lapack::gbtrf . | |
int64_t | lapack::gbcon (lapack::Norm norm, int64_t n, int64_t kl, int64_t ku, std::complex< float > const *AB, int64_t ldab, int64_t const *ipiv, float anorm, float *rcond) |
int64_t | lapack::gbequ (int64_t m, int64_t n, int64_t kl, int64_t ku, double const *AB, int64_t ldab, double *R, double *C, double *rowcnd, double *colcnd, double *amax) |
int64_t | lapack::gbequ (int64_t m, int64_t n, int64_t kl, int64_t ku, float const *AB, int64_t ldab, float *R, float *C, float *rowcnd, float *colcnd, float *amax) |
int64_t | lapack::gbequ (int64_t m, int64_t n, int64_t kl, int64_t ku, std::complex< double > const *AB, int64_t ldab, double *R, double *C, double *rowcnd, double *colcnd, double *amax) |
Computes row and column scalings intended to equilibrate an m-by-n band matrix A and reduce its condition number. | |
int64_t | lapack::gbequ (int64_t m, int64_t n, int64_t kl, int64_t ku, std::complex< float > const *AB, int64_t ldab, float *R, float *C, float *rowcnd, float *colcnd, float *amax) |
int64_t | lapack::gbequb (int64_t m, int64_t n, int64_t kl, int64_t ku, double const *AB, int64_t ldab, double *R, double *C, double *rowcnd, double *colcnd, double *amax) |
int64_t | lapack::gbequb (int64_t m, int64_t n, int64_t kl, int64_t ku, float const *AB, int64_t ldab, float *R, float *C, float *rowcnd, float *colcnd, float *amax) |
int64_t | lapack::gbequb (int64_t m, int64_t n, int64_t kl, int64_t ku, std::complex< double > const *AB, int64_t ldab, 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::gbequb (int64_t m, int64_t n, int64_t kl, int64_t ku, std::complex< float > const *AB, int64_t ldab, float *R, float *C, float *rowcnd, float *colcnd, float *amax) |
int64_t | lapack::gbrfs (lapack::Op trans, int64_t n, int64_t kl, int64_t ku, int64_t nrhs, double const *AB, int64_t ldab, double const *AFB, int64_t ldafb, int64_t const *ipiv, double const *B, int64_t ldb, double *X, int64_t ldx, double *ferr, double *berr) |
int64_t | lapack::gbrfs (lapack::Op trans, int64_t n, int64_t kl, int64_t ku, int64_t nrhs, float const *AB, int64_t ldab, float const *AFB, int64_t ldafb, int64_t const *ipiv, float const *B, int64_t ldb, float *X, int64_t ldx, float *ferr, float *berr) |
int64_t | lapack::gbrfs (lapack::Op trans, int64_t n, int64_t kl, int64_t ku, int64_t nrhs, std::complex< double > const *AB, int64_t ldab, std::complex< double > const *AFB, int64_t ldafb, 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 banded, and provides error bounds and backward error estimates for the solution. | |
int64_t | lapack::gbrfs (lapack::Op trans, int64_t n, int64_t kl, int64_t ku, int64_t nrhs, std::complex< float > const *AB, int64_t ldab, std::complex< float > const *AFB, int64_t ldafb, 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::gbtrf (int64_t m, int64_t n, int64_t kl, int64_t ku, double *AB, int64_t ldab, int64_t *ipiv) |
int64_t | lapack::gbtrf (int64_t m, int64_t n, int64_t kl, int64_t ku, float *AB, int64_t ldab, int64_t *ipiv) |
int64_t | lapack::gbtrf (int64_t m, int64_t n, int64_t kl, int64_t ku, std::complex< double > *AB, int64_t ldab, int64_t *ipiv) |
Computes an LU factorization of an m-by-n band matrix A using partial pivoting with row interchanges. | |
int64_t | lapack::gbtrf (int64_t m, int64_t n, int64_t kl, int64_t ku, std::complex< float > *AB, int64_t ldab, int64_t *ipiv) |
int64_t | lapack::gbtrs (lapack::Op trans, int64_t n, int64_t kl, int64_t ku, int64_t nrhs, double const *AB, int64_t ldab, int64_t const *ipiv, double *B, int64_t ldb) |
int64_t | lapack::gbtrs (lapack::Op trans, int64_t n, int64_t kl, int64_t ku, int64_t nrhs, float const *AB, int64_t ldab, int64_t const *ipiv, float *B, int64_t ldb) |
int64_t | lapack::gbtrs (lapack::Op trans, int64_t n, int64_t kl, int64_t ku, int64_t nrhs, std::complex< double > const *AB, int64_t ldab, int64_t const *ipiv, std::complex< double > *B, int64_t ldb) |
Solves a system of linear equations. | |
int64_t | lapack::gbtrs (lapack::Op trans, int64_t n, int64_t kl, int64_t ku, int64_t nrhs, std::complex< float > const *AB, int64_t ldab, int64_t const *ipiv, std::complex< float > *B, int64_t ldb) |
int64_t lapack::gbcon | ( | lapack::Norm | norm, |
int64_t | n, | ||
int64_t | kl, | ||
int64_t | ku, | ||
std::complex< double > const * | AB, | ||
int64_t | ldab, | ||
int64_t const * | ipiv, | ||
double | anorm, | ||
double * | rcond | ||
) |
Estimates the reciprocal of the condition number of a general band matrix A, in either the 1-norm or the infinity-norm, using the LU factorization computed by lapack::gbtrf
.
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>
.
[in] | norm | 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] | kl | The number of subdiagonals within the band of A. kl >= 0. |
[in] | ku | The number of superdiagonals within the band of A. ku >= 0. |
[in] | AB | The n-by-n band matrix AB, stored in an ldab-by-n array. Details of the LU factorization of the band matrix A, as computed by lapack::gbtrf . U is stored as an upper triangular band matrix with kl+ku superdiagonals in rows 1 to kl+ku+1, and the multipliers used during the factorization are stored in rows kl+ku+2 to 2*kl+ku+1. |
[in] | ldab | The leading dimension of the array AB. ldab >= 2*kl+ku+1. |
[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). |
[in] | anorm |
|
[out] | rcond | The reciprocal of the condition number of the matrix A, computed as rcond = 1/(norm(A) * norm(inv(A))). |
int64_t lapack::gbequ | ( | int64_t | m, |
int64_t | n, | ||
int64_t | kl, | ||
int64_t | ku, | ||
std::complex< double > const * | AB, | ||
int64_t | ldab, | ||
double * | R, | ||
double * | C, | ||
double * | rowcnd, | ||
double * | colcnd, | ||
double * | amax | ||
) |
Computes row and column scalings intended to equilibrate an m-by-n band 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>
.
[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] | kl | The number of subdiagonals within the band of A. kl >= 0. |
[in] | ku | The number of superdiagonals within the band of A. ku >= 0. |
[in] | AB | The n-by-n band matrix AB, stored in an ldab-by-n array. The band matrix A, stored in rows 1 to kl+ku+1. The j-th column of A is stored in the j-th column of the array AB as follows: AB(ku+1+i-j,j) = A(i,j) for max(1,j-ku) <= i <= min(m,j+kl). |
[in] | ldab | The leading dimension of the array AB. ldab >= kl+ku+1. |
[out] | R | The vector R of length m. If successful or return value > m, R contains the row scale factors for A. |
[out] | C | The vector C of length n. If successful, C contains the column scale factors for A. |
[out] | rowcnd | If 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] | colcnd | If 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] | amax | Absolute value of largest matrix element. If amax is very close to overflow or very close to underflow, the matrix should be scaled. |
int64_t lapack::gbequb | ( | int64_t | m, |
int64_t | n, | ||
int64_t | kl, | ||
int64_t | ku, | ||
std::complex< double > const * | AB, | ||
int64_t | ldab, | ||
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>
.
[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] | kl | The number of subdiagonals within the band of A. kl >= 0. |
[in] | ku | The number of superdiagonals within the band of A. ku >= 0. |
[in] | AB | The m-by-n band matrix AB, stored in an ldab-by-n array. On entry, the matrix A in band storage, in rows 1 to kl+ku+1. The j-th column of A is stored in the j-th column of the array AB as follows: AB(ku+1+i-j,j) = A(i,j) for max(1,j-ku) <= i <= min(n,j+kl) |
[in] | ldab | The leading dimension of the array A. ldab >= max(1,m). |
[out] | R | The vector R of length m. If successful or return value > m, R contains the row scale factors for A. |
[out] | C | The vector C of length n. If successful, C contains the column scale factors for A. |
[out] | rowcnd | If 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] | colcnd | If 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] | amax | Absolute value of largest matrix element. If amax is very close to overflow or very close to underflow, the matrix should be scaled. |
int64_t lapack::gbrfs | ( | lapack::Op | trans, |
int64_t | n, | ||
int64_t | kl, | ||
int64_t | ku, | ||
int64_t | nrhs, | ||
std::complex< double > const * | AB, | ||
int64_t | ldab, | ||
std::complex< double > const * | AFB, | ||
int64_t | ldafb, | ||
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 banded, 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 | The form of the system of equations:
|
[in] | n | The order of the matrix A. n >= 0. |
[in] | kl | The number of subdiagonals within the band of A. kl >= 0. |
[in] | ku | The number of superdiagonals within the band of A. ku >= 0. |
[in] | nrhs | The number of right hand sides, i.e., the number of columns of the matrices B and X. nrhs >= 0. |
[in] | AB | The n-by-n band matrix AB, stored in an ldab-by-n array. The original band matrix A, stored in rows 1 to kl+ku+1. The j-th column of A is stored in the j-th column of the array AB as follows: AB(ku+1+i-j,j) = A(i,j) for max(1,j-ku) <= i <= min(n,j+kl). |
[in] | ldab | The leading dimension of the array AB. ldab >= kl+ku+1. |
[in] | AFB | The n-by-n band matrix AFB, stored in an ldafb-by-n array. Details of the LU factorization of the band matrix A, as computed by lapack::gbtrf . U is stored as an upper triangular band matrix with kl+ku superdiagonals in rows 1 to kl+ku+1, and the multipliers used during the factorization are stored in rows kl+ku+2 to 2*kl+ku+1. |
[in] | ldafb | The leading dimension of the array AFB. ldafb >= 2*kl*ku+1. |
[in] | ipiv | The vector ipiv of length n. The pivot indices from lapack::gbtrf ; for 1 <= i <= n, row i of the matrix was interchanged with row ipiv(i). |
[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::gbtrs . 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::gbtrf | ( | int64_t | m, |
int64_t | n, | ||
int64_t | kl, | ||
int64_t | ku, | ||
std::complex< double > * | AB, | ||
int64_t | ldab, | ||
int64_t * | ipiv | ||
) |
Computes an LU factorization of an m-by-n band matrix A using partial pivoting with row interchanges.
This is the blocked version of the algorithm, calling Level 3 BLAS.
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] | kl | The number of subdiagonals within the band of A. kl >= 0. |
[in] | ku | The number of superdiagonals within the band of A. ku >= 0. |
[in,out] | AB | The n-by-n band matrix AB, stored in an ldab-by-n array. On entry, the matrix A in band storage, in rows kl+1 to 2*kl+ku+1; rows 1 to kl of the array need not be set. The j-th column of A is stored in the j-th column of the array AB as follows: AB(kl+ku+1+i-j,j) = A(i,j) for max(1,j-ku) <= i <= min(m,j+kl) On exit, details of the factorization: U is stored as an upper triangular band matrix with kl+ku superdiagonals in rows 1 to kl+ku+1, and the multipliers used during the factorization are stored in rows kl+ku+2 to 2*kl+ku+1. See below for further details. |
[in] | ldab | The leading dimension of the array AB. ldab >= 2*kl+ku+1. |
[out] | ipiv | The 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). |
The band storage scheme is illustrated by the following example, when m = n = 6, kl = 2, ku = 1:
On entry: On exit: * * * + + + * * * u14 u25 u36 * * + + + + * * u13 u24 u35 u46 * a12 a23 a34 a45 a56 * u12 u23 u34 u45 u56 a11 a22 a33 a44 a55 a66 u11 u22 u33 u44 u55 u66 a21 a32 a43 a54 a65 * m21 m32 m43 m54 m65 * a31 a42 a53 a64 * * m31 m42 m53 m64 * *
Array elements marked * are not used by the routine; elements marked
int64_t lapack::gbtrs | ( | lapack::Op | trans, |
int64_t | n, | ||
int64_t | kl, | ||
int64_t | ku, | ||
int64_t | nrhs, | ||
std::complex< double > const * | AB, | ||
int64_t | ldab, | ||
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 band matrix A using the LU factorization computed by lapack::gbtrf
.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[in] | trans | The form of the system of equations.
|
[in] | n | The order of the matrix A. n >= 0. |
[in] | kl | The number of subdiagonals within the band of A. kl >= 0. |
[in] | ku | The number of superdiagonals within the band of A. ku >= 0. |
[in] | nrhs | The number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0. |
[in] | AB | The n-by-n band matrix AB, stored in an ldab-by-n array. Details of the LU factorization of the band matrix A, as computed by lapack::gbtrf . U is stored as an upper triangular band matrix with kl+ku superdiagonals in rows 1 to kl+ku+1, and the multipliers used during the factorization are stored in rows kl+ku+2 to 2*kl+ku+1. |
[in] | ldab | The leading dimension of the array AB. ldab >= 2*kl+ku+1. |
[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). |
[in,out] | B | The 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] | ldb | The leading dimension of the array B. ldb >= max(1,n). |