LAPACK++ 2024.05.31
LAPACK C++ API
Loading...
Searching...
No Matches
Hermitian indefinite: Aasen's

Functions

int64_t lapack::hetrf_aa (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 Aasen's algorithm.
 
int64_t lapack::hetrf_aa (lapack::Uplo uplo, int64_t n, std::complex< float > *A, int64_t lda, int64_t *ipiv)
 
int64_t lapack::hetrs_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 hermitian matrix A using the factorization \(A = U T U^H\) or \(A = L T L^T\) computed by lapack::hetrf_aa.
 
int64_t lapack::hetrs_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)
 

Detailed Description

Function Documentation

◆ hetrf_aa()

int64_t lapack::hetrf_aa ( 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 Aasen's algorithm.

The form of the factorization is

\[ A = U T U^H \]

or

\[ A = L T L^H \]

where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and T is a hermitian 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, this is an alias for lapack::sytrf_aa. For complex symmetric matrices, see lapack::sytrf_aa.

Parameters
[in]uplo
  • lapack::Uplo::Upper: Upper triangle of A is stored;
  • lapack::Uplo::Lower: Lower triangle of A is stored.
[in]nThe order of the matrix A. n >= 0.
[in,out]AThe n-by-n matrix A, stored in an lda-by-n array.
  • On entry, the hermitian matrix A. If uplo = Upper, the leading n-by-n upper triangular part of A contains the upper triangular part of the matrix A, and the strictly lower triangular part of A is not referenced. If uplo = Lower, the leading n-by-n lower triangular part of A contains the lower triangular part of the matrix A, and the strictly upper triangular part of A is not referenced.
  • On exit, the tridiagonal matrix is stored in the diagonals and the subdiagonals of A just below (or above) the diagonals, and L is stored below (or above) the subdiaonals, when uplo is 'L' (or 'U').
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[out]ipivThe 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).
Returns
= 0: successful exit

◆ hetrs_aa()

int64_t lapack::hetrs_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 hermitian matrix A using the factorization \(A = U T U^H\) or \(A = L T L^T\) computed by lapack::hetrf_aa.

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>. For real matrices, this is an alias for lapack::sytrs_aa. For complex symmetric matrices, see lapack::sytrs_aa.

Parameters
[in]uploWhether the details of the factorization are stored as an upper or lower triangular matrix.
  • lapack::Uplo::Upper: Upper triangular, form is \(A = U T U^H;\)
  • lapack::Uplo::Lower: Lower triangular, form is \(A = L T L^H.\)
[in]nThe order of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0.
[in]AThe n-by-n matrix A, stored in an lda-by-n array. Details of factors computed by lapack::hetrf_aa.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in]ipivThe vector ipiv of length n. Details of the interchanges as computed by lapack::hetrf_aa.
[in,out]BThe 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]ldbThe leading dimension of the array B. ldb >= max(1,n).
Returns
= 0: successful exit