LAPACK++ 2024.05.31
LAPACK C++ API
|
Macros | |
#define | hetrf_3 hetrf_rk |
#define | hetri_3 hetri_rk |
#define | hetrs_3 hetrs_rk |
Functions | |
int64_t | lapack::hecon_rk (lapack::Uplo uplo, int64_t n, std::complex< double > const *A, int64_t lda, std::complex< double > const *E, 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 computed by lapack::hetrf_rk : | |
int64_t | lapack::hecon_rk (lapack::Uplo uplo, int64_t n, std::complex< float > const *A, int64_t lda, std::complex< float > const *E, int64_t const *ipiv, float anorm, float *rcond) |
int64_t | lapack::hetrf_rk (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, std::complex< double > *E, int64_t *ipiv) |
Computes the factorization of a Hermitian matrix A using the bounded Bunch-Kaufman (rook) diagonal pivoting method: | |
int64_t | lapack::hetrf_rk (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, std::complex< float > *E, int64_t *ipiv) |
int64_t | lapack::hetrf_rook (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, int64_t *ipiv) |
int64_t | lapack::hetrf_rook (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, int64_t *ipiv) |
int64_t | lapack::hetri_rk (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, std::complex< double > const *E, int64_t const *ipiv) |
Computes the inverse of a Hermitian indefinite matrix A using the factorization computed by lapack::hetrf_rk : | |
int64_t | lapack::hetri_rk (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, std::complex< float > const *E, int64_t const *ipiv) |
int64_t | lapack::hetrs_rk (lapack::Uplo uplo, int64_t n, int64_t nrhs, std::complex< double > const *A, int64_t lda, std::complex< double > const *E, 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 computed by lapack::hetrf_rk : | |
int64_t | lapack::hetrs_rk (lapack::Uplo uplo, int64_t n, int64_t nrhs, std::complex< float > const *A, int64_t lda, std::complex< float > const *E, int64_t const *ipiv, std::complex< float > *B, int64_t ldb) |
int64_t | lapack::hetrs_rook (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) |
int64_t | lapack::hetrs_rook (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) |
#define hetrf_3 hetrf_rk |
#define hetri_3 hetri_rk |
#define hetrs_3 hetrs_rk |
int64_t lapack::hecon_rk | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > const * | A, | ||
int64_t | lda, | ||
std::complex< double > const * | E, | ||
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 computed by lapack::hetrf_rk
:
\[ A = P U D U^H P^T \]
or
\[ A = P L D L^H P^T, \]
where U (or L) is unit upper (or lower) triangular matrix, \(U^H\) (or \(L^H\)) is the conjugate of U (or L), P is a permutation matrix, \(P^T\) is the transpose of P, and D is Hermitian and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.
An estimate is obtained for \(|| A^{-1} ||_1,\) and the reciprocal of the condition number is computed as \(\text{rcond} = 1 / (|| A ||_1 * || A^{-1} ||_1).\) This routine uses the BLAS-3 solver lapack::hetrs_rk
.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this is an alias for lapack::sycon_rk
. For complex symmetric matrices, see lapack::sycon_rk
.
[in] | uplo | Specifies 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. Diagonal of the block diagonal matrix D and factors U or L as computed by lapack::hetrf_rk :
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | E | The vector E of length n. On entry, contains the superdiagonal (or subdiagonal) elements of the Hermitian block diagonal matrix D with 1-by-1 or 2-by-2 diagonal blocks, where
|
[in] | ipiv | The vector ipiv of length n. Details of the interchanges and the block structure of D as determined by lapack::hetrf_rk . |
[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::hetrf_rk | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
std::complex< double > * | E, | ||
int64_t * | ipiv | ||
) |
Computes the factorization of a Hermitian matrix A using the bounded Bunch-Kaufman (rook) diagonal pivoting method:
\[ A = P U D U^H P^T \]
or
\[ A = P L D L^H P^T, \]
where U (or L) is unit upper (or lower) triangular matrix, \(U^H\) (or \(L^H\)) is the conjugate of U (or L), P is a permutation matrix, \(P^T\) is the transpose of P, 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. For more information see Further Details section.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this is an alias for lapack::sytrf_rk
. For complex symmetric matrices, see lapack::sytrf_rk
.
lapack::hetrf_rook
.[in] | uplo | Whether the upper or lower triangular part of the Hermitian matrix A is stored:
|
[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.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[out] | E | The vector E of length n. On exit, contains the superdiagonal (or subdiagonal) elements of the Hermitian block diagonal matrix D with 1-by-1 or 2-by-2 diagonal blocks, where
|
[out] | ipiv | The vector ipiv of length n. ipiv describes the permutation matrix P in the factorization of matrix A as follows. The absolute value of ipiv(k) represents the index of row and column that were interchanged with the k-th row and column. The value of uplo describes the order in which the interchanges were applied. Also, the sign of ipiv represents the block structure of the Hermitian block diagonal matrix D with 1-by-1 or 2-by-2 diagonal blocks which correspond to 1 or 2 interchanges at each factorization step. For more info see Further Details section.
|
Therefore D(i,i) is exactly zero, and superdiagonal elements of column i of U (or subdiagonal elements of column i of L ) are all zeros. The factorization has been completed, but the block diagonal matrix D is exactly singular, and division by zero will occur if it is used to solve a system of equations.
NOTE: info only stores the first occurrence of a singularity, any subsequent occurrence of singularity is not stored in info even though the factorization always completes.
int64_t lapack::hetrf_rook | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
int64_t * | ipiv | ||
) |
int64_t lapack::hetri_rk | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
std::complex< double > const * | E, | ||
int64_t const * | ipiv | ||
) |
Computes the inverse of a Hermitian indefinite matrix A using the factorization computed by lapack::hetrf_rk
:
\[ A = P U D U^H P^T \]
or
\[ A = P L D L^H P^T, \]
where U (or L) is unit upper (or lower) triangular matrix, \(U^H\) (or \(L^H\)) is the conjugate of U (or L), P is a permutation matrix, \(P^T\) is the transpose of P, 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 is an alias for lapack::sytri_rk
. For complex symmetric matrices, see lapack::sytri_rk
.
Note: LAPACK++ uses the name hetri_rk
instead of LAPACK's hetri_3
, for consistency with hesv_rk
, hetrf_rk
, etc.
lapack::hetri_rook
.[in] | uplo | Specifies 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.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | E | The vector E of length n. On entry, contains the superdiagonal (or subdiagonal) elements of the Hermitian block diagonal matrix D with 1-by-1 or 2-by-2 diagonal blocks, where
|
[in] | ipiv | The vector ipiv of length n. Details of the interchanges and the block structure of D as determined by lapack::hetrf_rk . |
int64_t lapack::hetrs_rk | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
int64_t | nrhs, | ||
std::complex< double > const * | A, | ||
int64_t | lda, | ||
std::complex< double > const * | E, | ||
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 computed by lapack::hetrf_rk
:
\[ A = P U D U^H P^T \]
or
\[ A = P L D L^H P^T, \]
where U (or L) is unit upper (or lower) triangular matrix, \(U^H\) (or \(L^H\)) is the conjugate of U (or L), P is a permutation matrix, \(P^T\) is the transpose of P, and D is Hermitian and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.
This algorithm is using Level 3 BLAS.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, this is an alias for lapack::sytrs_rk
. For complex symmetric matrices, see lapack::sytrs_rk
.
Note: LAPACK++ uses the name hetrs_rk
instead of LAPACK's hetrs_3
, for consistency with hesv_rk
, hetrf_rk
, etc.
lapack::hetrf_rook
.[in] | uplo | Specifies 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. Diagonal of the block diagonal matrix D and factors U or L as computed by lapack::hetrf_rk :
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[in] | E | The vector E of length n. On entry, contains the superdiagonal (or subdiagonal) elements of the Hermitian block diagonal matrix D with 1-by-1 or 2-by-2 diagonal blocks, where
|
[in] | ipiv | The vector ipiv of length n. Details of the interchanges and the block structure of D as determined by lapack::hetrf_rk . |
[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::hetrs_rook | ( | 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 | ||
) |