LAPACK++ 2024.05.31
LAPACK C++ API
|
Macros | |
#define | sytrf_3 sytrf_rk |
#define | sytri_3 sytri_rk |
#define | sytrs_3 sytrs_rk |
Functions | |
int64_t | lapack::sycon_rk (lapack::Uplo uplo, int64_t n, double const *A, int64_t lda, double const *E, int64_t const *ipiv, double anorm, double *rcond) |
int64_t | lapack::sycon_rk (lapack::Uplo uplo, int64_t n, float const *A, int64_t lda, float const *E, int64_t const *ipiv, float anorm, float *rcond) |
int64_t | lapack::sycon_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 symmetric matrix A using the factorization computed by lapack::sytrf_rk : | |
int64_t | lapack::sycon_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::sytrf_rk (lapack::Uplo uplo, int64_t n, double *A, int64_t lda, double *E, int64_t *ipiv) |
int64_t | lapack::sytrf_rk (lapack::Uplo uplo, int64_t n, float *A, int64_t lda, float *E, int64_t *ipiv) |
int64_t | lapack::sytrf_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 symmetric matrix A using the bounded Bunch-Kaufman (rook) diagonal pivoting method: | |
int64_t | lapack::sytrf_rk (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, std::complex< float > *E, int64_t *ipiv) |
int64_t | lapack::sytrf_rook (lapack::Uplo uplo, int64_t n, double *A, int64_t lda, int64_t *ipiv) |
int64_t | lapack::sytrf_rook (lapack::Uplo uplo, int64_t n, float *A, int64_t lda, int64_t *ipiv) |
int64_t | lapack::sytrf_rook (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, int64_t *ipiv) |
int64_t | lapack::sytrf_rook (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, int64_t *ipiv) |
int64_t | lapack::sytri_rk (lapack::Uplo uplo, int64_t n, double *A, int64_t lda, double const *E, int64_t const *ipiv) |
int64_t | lapack::sytri_rk (lapack::Uplo uplo, int64_t n, float *A, int64_t lda, float const *E, int64_t const *ipiv) |
int64_t | lapack::sytri_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 symmetric indefinite matrix A using the factorization computed by lapack::sytrf_rk : | |
int64_t | lapack::sytri_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::sytrs_rk (lapack::Uplo uplo, int64_t n, int64_t nrhs, double const *A, int64_t lda, double const *E, int64_t const *ipiv, double *B, int64_t ldb) |
int64_t | lapack::sytrs_rk (lapack::Uplo uplo, int64_t n, int64_t nrhs, float const *A, int64_t lda, float const *E, int64_t const *ipiv, float *B, int64_t ldb) |
int64_t | lapack::sytrs_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 symmetric matrix A using the factorization computed by lapack::sytrf_rk : | |
int64_t | lapack::sytrs_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::sytrs_rook (lapack::Uplo uplo, int64_t n, int64_t nrhs, double const *A, int64_t lda, int64_t const *ipiv, double *B, int64_t ldb) |
int64_t | lapack::sytrs_rook (lapack::Uplo uplo, int64_t n, int64_t nrhs, float const *A, int64_t lda, int64_t const *ipiv, float *B, int64_t ldb) |
int64_t | lapack::sytrs_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::sytrs_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 sytrf_3 sytrf_rk |
#define sytri_3 sytri_rk |
#define sytrs_3 sytrs_rk |
int64_t lapack::sycon_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 symmetric matrix A using the factorization computed by lapack::sytrf_rk
:
\[ A = P U D U^T P^T, \]
or
\[ A = P L D L^T P^T, \]
where U (or L) is unit upper (or lower) triangular matrix, \(U^T\) (or \(L^T\)) is the transpose of U (or L), P is a permutation matrix, \(P^T\) is the transpose of P, and D is symmetric 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, lapack::hecon_rk
is an alias for this. For complex Hermitian matrices, see lapack::hecon_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::sytrf_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 symmetric 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::sytrf_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::sytrf_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 symmetric matrix A using the bounded Bunch-Kaufman (rook) diagonal pivoting method:
\[ A = P U D U^T P^T \]
or
\[ A = P L D L^T P^T, \]
where U (or L) is unit upper (or lower) triangular matrix, \(U^T\) (or \(L^T\)) is the transpose of U (or L), P is a permutation matrix, \(P^T\) is the transpose of P, and D is symmetric 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, lapack::hetrf_rk
is an alias for this. For complex Hermitian matrices, see lapack::hetrf_rk
.
lapack::sytrf_rook
.[in] | uplo | Whether the upper or lower triangular part of the symmetric 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 symmetric 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 symmetric 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::sytrf_rook | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
int64_t * | ipiv | ||
) |
int64_t lapack::sytri_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 symmetric indefinite matrix A using the factorization computed by lapack::sytrf_rk
:
\[ A = P U D U^T P^T \]
or
\[ A = P L D L^T P^T, \]
where U (or L) is unit upper (or lower) triangular matrix, \(U^T\) (or \(L^T\)) is the transpose of U (or L), P is a permutation matrix, \(P^T\) is the transpose of P, and D is symmetric 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, lapack::hetri_rk
is an alias for this. For complex Hermitian matrices, see lapack::hetri_rk
.
Note: LAPACK++ uses the name sytri_rk
instead of LAPACK's sytri_3
, for consistency with sysv_rk
, sytrf_rk
, etc.
lapack::sytri_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 symmetric 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::sytrf_rk . |
int64_t lapack::sytrs_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 symmetric matrix A using the factorization computed by lapack::sytrf_rk
:
\[ A = P U D (U^T) (P^T)\f$ or \f$A = P L D (L^T) (P^T), \]
where U (or L) is unit upper (or lower) triangular matrix, \(U^T\) (or \(L^T\)) is the transpose of U (or L), P is a permutation matrix, \(P^T\) is the transpose of P, and D is symmetric 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, lapack::hetrs_rk
is an alias for this. For complex Hermitian matrices, see lapack::hetrs_rk
.
Note: LAPACK++ uses the name sytrs_rk
instead of LAPACK's sytrs_3
, for consistency with sysv_rk
, sytrf_rk
, etc.
lapack::sytrs_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::sytrf_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 symmetric 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::sytrf_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::sytrs_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 | ||
) |