LAPACK++ 2024.05.31
LAPACK C++ API
|
Functions | |
int64_t | lapack::hecon (lapack::Uplo uplo, int64_t n, std::complex< double > const *A, int64_t lda, int64_t const *ipiv, double anorm, double *rcond) |
Estimates the reciprocal of the condition number (in the 1-norm) of a Hermitian matrix A using the factorization \(A = U D U^H\) or \(A = L D L^H\) computed by lapack::hetrf . | |
int64_t | lapack::hecon (lapack::Uplo uplo, int64_t n, std::complex< float > const *A, int64_t lda, int64_t const *ipiv, float anorm, float *rcond) |
int64_t | lapack::heequb (lapack::Uplo uplo, int64_t n, std::complex< double > const *A, int64_t lda, double *S, double *scond, double *amax) |
Computes row and column scalings intended to equilibrate a Hermitian matrix A (with respect to the Euclidean norm) and reduce its condition number. | |
int64_t | lapack::heequb (lapack::Uplo uplo, int64_t n, std::complex< float > const *A, int64_t lda, float *S, float *scond, float *amax) |
int64_t | lapack::herfs (lapack::Uplo uplo, 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 when the coefficient matrix is Hermitian indefinite, and provides error bounds and backward error estimates for the solution. | |
int64_t | lapack::herfs (lapack::Uplo uplo, 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) |
void | lapack::heswapr (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, int64_t i1, int64_t i2) |
Applies an elementary permutation on the rows and the columns of a Hermitian matrix. | |
void | lapack::heswapr (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, int64_t i1, int64_t i2) |
int64_t | lapack::hetrf (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, int64_t *ipiv) |
Computes the factorization of a Hermitian matrix A using the Bunch-Kaufman diagonal pivoting method. | |
int64_t | lapack::hetrf (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, int64_t *ipiv) |
int64_t | lapack::hetri (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, int64_t const *ipiv) |
Computes the inverse of a Hermitian indefinite matrix A using the factorization \(A = U D U^H\) or \(A = L D L^H\) computed by lapack::hetrf . | |
int64_t | lapack::hetri (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, int64_t const *ipiv) |
int64_t | lapack::hetri2 (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, int64_t const *ipiv) |
Computes the inverse of a Hermitian indefinite matrix A using the factorization \(A = U D U^H\) or \(A = L D L^H\) computed by lapack::hetrf . | |
int64_t | lapack::hetri2 (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, int64_t const *ipiv) |
int64_t | lapack::hetrs (lapack::Uplo uplo, 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\) with a Hermitian matrix A using the factorization \(A = U D U^H\) or \(A = L D L^H\) computed by lapack::hetrf . | |
int64_t | lapack::hetrs (lapack::Uplo uplo, 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) |
int64_t | lapack::hetrs2 (lapack::Uplo uplo, 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\) with a Hermitian matrix A using the factorization \(A = U D U^H\) or \(A = L D L^H\) computed by lapack::hetrf and converted by lapack::syconv . | |
int64_t | lapack::hetrs2 (lapack::Uplo uplo, 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) |
int64_t lapack::hecon | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > const * | A, | ||
int64_t | lda, | ||
int64_t const * | ipiv, | ||
double | anorm, | ||
double * | rcond | ||
) |
Estimates the reciprocal of the condition number (in the 1-norm) of a Hermitian matrix A using the factorization \(A = U D U^H\) or \(A = L D L^H\) computed by lapack::hetrf
.
An estimate is obtained for \(|| A^{-1} ||_1,\) and the reciprocal of the condition number is computed as \(\text{rcond} = 1 / (||A||_1 \cdot || A^{-1} ||_1).\)
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this in an alias for lapack::sycon
. For complex symmetric matrices, see lapack::sycon
.
[in] | uplo | Whether the details of the factorization are stored as an upper or lower triangular matrix.
|
[in] | n | The order of the matrix A. n >= 0. |
[in] | A | The n-by-n matrix A, stored in an lda-by-n array. The block diagonal matrix D and the multipliers used to obtain the factor U or L as computed by lapack::hetrf . |
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | ipiv | The vector ipiv of length n. Details of the interchanges and the block structure of D as determined by lapack::hetrf . |
[in] | anorm | The 1-norm of the original matrix A. |
[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 of \(A^{-1}\) computed in this routine. |
int64_t lapack::heequb | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > const * | A, | ||
int64_t | lda, | ||
double * | S, | ||
double * | scond, | ||
double * | amax | ||
) |
Computes row and column scalings intended to equilibrate a Hermitian matrix A (with respect to the Euclidean norm) and reduce its condition number.
The scale factors S are computed by the BIN algorithm (see references) so that the scaled matrix B with elements \(B_{i,j} = S_{i} A_{i,j} S_{j}\) has a condition number within a factor n of the smallest possible condition number over all possible diagonal scalings.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this in an alias for lapack::syequb
. For complex symmetric matrices, see lapack::syequb
.
[in] | uplo |
|
[in] | n | The order of the matrix A. n >= 0. |
[in] | A | The n-by-n matrix A, stored in an lda-by-n array. The n-by-n Hermitian matrix whose scaling factors are to be computed. |
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[out] | S | The vector S of length n. If successful, S contains the scale factors for A. |
[out] | scond | If successful, S contains the ratio of the smallest S(i) to the largest S(i). If scond >= 0.1 and amax is neither too large nor too small, it is not worth scaling by S. |
[out] | amax | Largest absolute value of any matrix element. If amax is very close to overflow or very close to underflow, the matrix should be scaled. |
int64_t lapack::herfs | ( | lapack::Uplo | uplo, |
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 when the coefficient matrix is Hermitian indefinite, 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>
. For real matrices, this in an alias for lapack::syrfs
. For complex symmetric matrices, see lapack::syrfs
.
[in] | uplo |
|
[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 matrices B and X. nrhs >= 0. |
[in] | A | The n-by-n matrix A, stored in an lda-by-n array. The Hermitian matrix A.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | AF | The n-by-n matrix AF, stored in an ldaf-by-n array. The factored form of the matrix A. AF contains the block diagonal matrix D and the multipliers used to obtain the factor U or L from the factorization \(A = U D U^H\) or \(A = L D L^H\) as computed by lapack::hetrf . |
[in] | ldaf | The leading dimension of the array AF. ldaf >= max(1,n). |
[in] | ipiv | The vector ipiv of length n. Details of the interchanges and the block structure of D as determined by lapack::hetrf . |
[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::hetrs . 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). |
void lapack::heswapr | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
int64_t | i1, | ||
int64_t | i2 | ||
) |
Applies an elementary permutation on the rows and the columns of a Hermitian matrix.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this in an alias for lapack::syswapr
. For complex symmetric matrices, see lapack::syswapr
.
[in] | uplo | Whether the details of the factorization are stored as an upper or lower triangular matrix.
|
[in] | n | The order of the matrix A. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. TODO (the LAPACK documentation seems wrong). |
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | i1 | Index of the first row to swap |
[in] | i2 | Index of the second row to swap |
int64_t lapack::hetrf | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
int64_t * | ipiv | ||
) |
Computes the factorization of a Hermitian matrix A using the Bunch-Kaufman diagonal pivoting method.
The form of the factorization is
\[ A = U D U^H \]
or
\[ A = L D L^H \]
where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and D is Hermitian and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.
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>
. For real matrices, this in an alias for lapack::sytrf
. For complex symmetric matrices, see lapack::sytrf
.
[in] | uplo |
|
[in] | n | The order of the matrix A. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. On entry, the Hermitian matrix A.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[out] | ipiv | The vector ipiv of length n. Details of the interchanges and the block structure of D.
|
If uplo = Upper, then \(A = U D U^H,\) where
\[ U = P(n) U(n) \dots P(k) U(k) \dots, \]
i.e., U is a product of terms \(P(k) U(k),\) where k decreases from n to 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by ipiv(k), and U(k) is a unit upper triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I v 0 ) k-s U(k) = ( 0 I 0 ) s ( 0 0 I ) n-k k-s s n-k
If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k). If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k), and A(k,k), and v overwrites A(1:k-2,k-1:k).
If uplo = Lower, then \(A = L D L^H,\) where
\[ L = P(1) L(1) \dots P(k) L(k) \dots, \]
i.e., L is a product of terms \(P(k) L(k),\) where k increases from 1 to n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1 and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as defined by ipiv(k), and L(k) is a unit lower triangular matrix, such that if the diagonal block D(k) is of order s (s = 1 or 2), then
( I 0 0 ) k-1 L(k) = ( 0 I 0 ) s ( 0 v I ) n-k-s+1 k-1 s n-k-s+1
If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k). If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k), and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
int64_t lapack::hetri | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
int64_t const * | ipiv | ||
) |
Computes the inverse of a Hermitian indefinite matrix A using the factorization \(A = U D U^H\) or \(A = L D L^H\) computed by lapack::hetrf
.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this in an alias for lapack::sytri
. For complex symmetric matrices, see lapack::sytri
.
[in] | uplo | Whether the details of the factorization are stored as an upper or lower triangular matrix.
|
[in] | n | The order of the matrix A. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. On entry, the block diagonal matrix D and the multipliers used to obtain the factor U or L as computed by lapack::hetrf . On successful exit, the (Hermitian) inverse of the original matrix.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | ipiv | The vector ipiv of length n. Details of the interchanges and the block structure of D as determined by lapack::hetrf . |
int64_t lapack::hetri2 | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
int64_t const * | ipiv | ||
) |
Computes the inverse of a Hermitian indefinite matrix A using the factorization \(A = U D U^H\) or \(A = L D L^H\) computed by lapack::hetrf
.
hetri2
sets the leading dimension of the workspace before calling hetri2x
that actually computes the inverse.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this in an alias for lapack::sytri2
. For complex symmetric matrices, see lapack::sytri2
.
[in] | uplo | Whether the details of the factorization are stored as an upper or lower triangular matrix.
|
[in] | n | The order of the matrix A. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. On entry, the block diagonal matrix D and the multipliers used to obtain the factor U or L as computed by lapack::hetrf . On successful exit, the (Hermitian) inverse of the original matrix.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | ipiv | The vector ipiv of length n. Details of the interchanges and the block structure of D as determined by lapack::hetrf . |
int64_t lapack::hetrs | ( | lapack::Uplo | uplo, |
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\) with a Hermitian matrix A using the factorization \(A = U D U^H\) or \(A = L D L^H\) computed by lapack::hetrf
.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this in an alias for lapack::sytrs
. For complex symmetric matrices, see lapack::sytrs
.
[in] | uplo | Whether the details of the factorization are stored as an upper or lower triangular matrix.
|
[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] | A | The n-by-n matrix A, stored in an lda-by-n array. The block diagonal matrix D and the multipliers used to obtain the factor U or L as computed by lapack::hetrf . |
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | ipiv | The vector ipiv of length n. Details of the interchanges and the block structure of D as determined by lapack::hetrf . |
[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). |
int64_t lapack::hetrs2 | ( | lapack::Uplo | uplo, |
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\) with a Hermitian matrix A using the factorization \(A = U D U^H\) or \(A = L D L^H\) computed by lapack::hetrf
and converted by lapack::syconv
.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this in an alias for lapack::sytrs2
. For complex symmetric matrices, see lapack::sytrs2
.
[in] | uplo | Whether the details of the factorization are stored as an upper or lower triangular matrix.
|
[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] | A | The n-by-n matrix A, stored in an lda-by-n array. The block diagonal matrix D and the multipliers used to obtain the factor U or L as computed by lapack::hetrf . |
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | ipiv | The vector ipiv of length n. Details of the interchanges and the block structure of D as determined by lapack::hetrf . |
[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). |