LAPACK++ 2024.05.31
LAPACK C++ API
|
Functions | |
int64_t | lapack::lauum (lapack::Uplo uplo, int64_t n, double *A, int64_t lda) |
int64_t | lapack::lauum (lapack::Uplo uplo, int64_t n, float *A, int64_t lda) |
int64_t | lapack::lauum (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda) |
Computes the product \(U U^H\) or \(L^H L,\) where the triangular factor U or L is stored in the upper or lower triangular part of the array A. | |
int64_t | lapack::lauum (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda) |
int64_t | lapack::pocon (lapack::Uplo uplo, int64_t n, double const *A, int64_t lda, double anorm, double *rcond) |
int64_t | lapack::pocon (lapack::Uplo uplo, int64_t n, float const *A, int64_t lda, float anorm, float *rcond) |
int64_t | lapack::pocon (lapack::Uplo uplo, int64_t n, std::complex< double > const *A, int64_t lda, double anorm, double *rcond) |
Estimates the reciprocal of the condition number (in the 1-norm) of a Hermitian positive definite matrix using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::potrf . | |
int64_t | lapack::pocon (lapack::Uplo uplo, int64_t n, std::complex< float > const *A, int64_t lda, float anorm, float *rcond) |
int64_t | lapack::poequ (int64_t n, double const *A, int64_t lda, double *S, double *scond, double *amax) |
int64_t | lapack::poequ (int64_t n, float const *A, int64_t lda, float *S, float *scond, float *amax) |
int64_t | lapack::poequ (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 positive definite matrix A and reduce its condition number (with respect to the two-norm). | |
int64_t | lapack::poequ (int64_t n, std::complex< float > const *A, int64_t lda, float *S, float *scond, float *amax) |
int64_t | lapack::poequb (int64_t n, double const *A, int64_t lda, double *S, double *scond, double *amax) |
int64_t | lapack::poequb (int64_t n, float const *A, int64_t lda, float *S, float *scond, float *amax) |
int64_t | lapack::poequb (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 positive definite matrix A and reduce its condition number (with respect to the two-norm). | |
int64_t | lapack::poequb (int64_t n, std::complex< float > const *A, int64_t lda, float *S, float *scond, float *amax) |
int64_t | lapack::porfs (lapack::Uplo uplo, int64_t n, int64_t nrhs, double const *A, int64_t lda, double const *AF, int64_t ldaf, double const *B, int64_t ldb, double *X, int64_t ldx, double *ferr, double *berr) |
int64_t | lapack::porfs (lapack::Uplo uplo, int64_t n, int64_t nrhs, float const *A, int64_t lda, float const *AF, int64_t ldaf, float const *B, int64_t ldb, float *X, int64_t ldx, float *ferr, float *berr) |
int64_t | lapack::porfs (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, 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 positive definite, and provides error bounds and backward error estimates for the solution. | |
int64_t | lapack::porfs (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, std::complex< float > const *B, int64_t ldb, std::complex< float > *X, int64_t ldx, float *ferr, float *berr) |
int64_t | lapack::potf2 (lapack::Uplo uplo, int64_t n, double *A, int64_t lda) |
int64_t | lapack::potf2 (lapack::Uplo uplo, int64_t n, float *A, int64_t lda) |
int64_t | lapack::potf2 (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A. | |
int64_t | lapack::potf2 (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda) |
int64_t | lapack::potrf (lapack::Uplo uplo, int64_t n, double *A, int64_t lda) |
int64_t | lapack::potrf (lapack::Uplo uplo, int64_t n, float *A, int64_t lda) |
int64_t | lapack::potrf (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A. | |
int64_t | lapack::potrf (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda) |
int64_t | lapack::potrf2 (lapack::Uplo uplo, int64_t n, double *A, int64_t lda) |
int64_t | lapack::potrf2 (lapack::Uplo uplo, int64_t n, float *A, int64_t lda) |
int64_t | lapack::potrf2 (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A using the recursive algorithm. | |
int64_t | lapack::potrf2 (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda) |
int64_t | lapack::potri (lapack::Uplo uplo, int64_t n, double *A, int64_t lda) |
int64_t | lapack::potri (lapack::Uplo uplo, int64_t n, float *A, int64_t lda) |
int64_t | lapack::potri (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda) |
Computes the inverse of a Hermitian positive definite matrix A using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::potrf . | |
int64_t | lapack::potri (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda) |
int64_t | lapack::potrs (lapack::Uplo uplo, int64_t n, int64_t nrhs, double const *A, int64_t lda, double *B, int64_t ldb) |
int64_t | lapack::potrs (lapack::Uplo uplo, int64_t n, int64_t nrhs, float const *A, int64_t lda, float *B, int64_t ldb) |
int64_t | lapack::potrs (lapack::Uplo uplo, int64_t n, int64_t nrhs, std::complex< double > const *A, int64_t lda, std::complex< double > *B, int64_t ldb) |
Solves a system of linear equations \(A X = B\) with a Hermitian positive definite matrix A using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::potrf . | |
int64_t | lapack::potrs (lapack::Uplo uplo, int64_t n, int64_t nrhs, std::complex< float > const *A, int64_t lda, std::complex< float > *B, int64_t ldb) |
int64_t lapack::lauum | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda | ||
) |
Computes the product \(U U^H\) or \(L^H L,\) where the triangular factor U or L is stored in the upper or lower triangular part of the array A.
If uplo = Upper then the upper triangle of the result is stored, overwriting the factor U in A.
If uplo = Lower then the lower triangle of the result is stored, overwriting the factor L in A.
This is the blocked form of the algorithm, calling Level 3 BLAS.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[in] | uplo | Whether the triangular factor stored in the array A is upper or lower triangular:
|
[in] | n | The order of the triangular factor U or L. n >= 0. |
[in,out] | A | The n-by-n matrix A, stored in an lda-by-n array. On entry, the triangular factor U or L. On exit,
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
int64_t lapack::pocon | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > const * | A, | ||
int64_t | lda, | ||
double | anorm, | ||
double * | rcond | ||
) |
Estimates the reciprocal of the condition number (in the 1-norm) of a Hermitian positive definite matrix using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::potrf
.
An estimate is obtained for norm(inv(A)), and the reciprocal of the condition number is computed as rcond = 1 / (anorm * norm(inv(A))).
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[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 triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H\), as computed by lapack::potrf . |
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | anorm | The 1-norm (or infinity-norm) of the Hermitian 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 inv(A) computed in this routine. |
int64_t lapack::poequ | ( | 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 positive definite matrix A and reduce its condition number (with respect to the two-norm).
S contains the scale factors, \(S_{i} = 1/\sqrt{ A_{i,i} },\) chosen so that the scaled matrix B with elements \(B_{i,j} = S_{i} A_{i,j} S_{j}\) has ones on the diagonal. This choice of S puts the condition number of B 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>
.
[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 positive definite matrix whose scaling factors are to be computed. Only the diagonal elements of A are referenced. |
[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 | 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::poequb | ( | 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 positive definite matrix A and reduce its condition number (with respect to the two-norm).
S contains the scale factors, \(S_{i} = 1/\sqrt{ A_{i,i} },\) chosen so that the scaled matrix B with elements \(B_{i,j} = S_{i} A_{i,j} S_{j}\) has ones on the diagonal. This choice of S puts the condition number of B within a factor n of the smallest possible condition number over all possible diagonal scalings.
This routine differs from lapack::poequ
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 diagonal entries 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] | 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 positive definite matrix whose scaling factors are to be computed. Only the diagonal elements of A are referenced. |
[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 | 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::porfs | ( | 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, | ||
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 positive definite, 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] | 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 triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H\), as computed by lapack::potrf . |
[in] | ldaf | The leading dimension of the array AF. ldaf >= max(1,n). |
[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::potrs . 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::potf2 | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda | ||
) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A.
The factorization has the form \(A = U^H U,\) if uplo = Upper, or \(A = L L^H,\) if uplo = Lower, where \(U\) is an upper triangular matrix and \(L\) is a lower triangular matrix.
This is the unblocked version of the algorithm, calling Level 2 BLAS.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[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). |
int64_t lapack::potrf | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda | ||
) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A.
The factorization has the form \(A = U^H U,\) if uplo = Upper, or \(A = L L^H,\) if uplo = Lower, where \(U\) is an upper triangular matrix and \(L\) is a lower triangular matrix.
This is the block version of the algorithm, calling Level 3 BLAS.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[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). |
int64_t lapack::potrf2 | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda | ||
) |
Computes the Cholesky factorization of a Hermitian positive definite matrix A using the recursive algorithm.
The factorization has the form \(A = U^H U,\) if uplo = Upper, or \(A = L L^H,\) if uplo = Lower, where U is an upper triangular matrix and L is lower triangular.
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 = n/2 and n2 = n-n1. The subroutine calls itself to factor \(A_{11},\) updates and scales \(A_{21}\) or \(A_{12},\) updates \(A_{22},\) and calls itself to factor \(A_{22}.\)
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[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). |
int64_t lapack::potri | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda | ||
) |
Computes the inverse of a Hermitian positive definite matrix A using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::potrf
.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[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 triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H,\) as computed by lapack::potrf . On exit, the upper or lower triangle of the (Hermitian) inverse of A, overwriting the input factor U or L. |
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
int64_t lapack::potrs | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
int64_t | nrhs, | ||
std::complex< double > const * | A, | ||
int64_t | lda, | ||
std::complex< double > * | B, | ||
int64_t | ldb | ||
) |
Solves a system of linear equations \(A X = B\) with a Hermitian positive definite matrix A using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::potrf
.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
.
[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 matrix B. nrhs >= 0. |
[in] | A | The n-by-n matrix A, stored in an lda-by-n array. The triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H\), as computed by lapack::potrf . |
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[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). |