LAPACK++ 2024.05.31
LAPACK C++ API
|
Functions | |
int64_t | lapack::sytrf_aa (lapack::Uplo uplo, int64_t n, double *A, int64_t lda, int64_t *ipiv) |
int64_t | lapack::sytrf_aa (lapack::Uplo uplo, int64_t n, float *A, int64_t lda, int64_t *ipiv) |
int64_t | lapack::sytrf_aa (lapack::Uplo uplo, int64_t n, std::complex< double > *A, int64_t lda, int64_t *ipiv) |
Computes the factorization of a symmetric matrix A using the Aasen's algorithm. | |
int64_t | lapack::sytrf_aa (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, int64_t *ipiv) |
int64_t | lapack::sytrs_aa (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_aa (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_aa (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 symmetric matrix A using the factorization \(A = U T U^T\) or \(A = L T L^T\) computed by lapack::sytrf_aa . | |
int64_t | lapack::sytrs_aa (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::sytrf_aa | ( | lapack::Uplo | uplo, |
int64_t | n, | ||
std::complex< double > * | A, | ||
int64_t | lda, | ||
int64_t * | ipiv | ||
) |
Computes the factorization of a symmetric matrix A using the Aasen's algorithm.
The form of the factorization is
\[ A = U T U^T \]
or
\[ A = L T L^T \]
where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and T is a symmetric tridiagonal matrix.
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::hetrf_aa
is an alias for this. For complex Hermitian matrices, see lapack::hetrf_aa
.
[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.
|
[in] | lda | The leading dimension of the array A. lda >= max(1,n). |
[out] | ipiv | The vector ipiv of length n. On exit, it contains the details of the interchanges, i.e., the row and column k of A were interchanged with the row and column ipiv(k). |
int64_t lapack::sytrs_aa | ( | 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 symmetric matrix A using the factorization \(A = U T U^T\) or \(A = L T L^T\) computed by lapack::sytrf_aa
.
Overloaded versions are available for float
, double
, std::complex<float>
, and std::complex<double>
. For real matrices, lapack::hetrs_aa
is an alias for this. For complex Hermitian matrices, see lapack::hetrs_aa
.
[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. Details of factors computed by lapack::sytrf_aa . |
[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 as computed by lapack::sytrf_aa . |
[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). |