hesv_aa#

Functions

void chesv_aa(
    const char*          uplo,
    const INT            n,
    const INT            nrhs,
          c64*  restrict A,
    const INT            lda,
          INT*  restrict ipiv,
          c64*  restrict B,
    const INT            ldb,
          c64*  restrict work,
    const INT            lwork,
          INT*           info
);
void chesv_aa(const char *uplo, const INT n, const INT nrhs, c64 *restrict A, const INT lda, INT *restrict ipiv, c64 *restrict B, const INT ldb, c64 *restrict work, const INT lwork, INT *info)#

CHESV_AA computes the solution to a complex system of linear equations A * X = B, where A is an N-by-N Hermitian matrix and X and B are N-by-NRHS matrices.

Aasen’s algorithm is used to factor A as A = U**H * T * U, if UPLO = ‘U’, or A = L * T * L**H, if UPLO = ‘L’, where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and T is Hermitian and tridiagonal. The factored form of A is then used to solve the system of equations A * X = B.

On exit, if info = 0, the tridiagonal matrix T and the multipliers used to obtain the factor U or L from the factorization A = U**H*T*U or A = L*T*L**H as computed by CHETRF_AA.

Parameters

in
uplo

= ‘U’: Upper triangle of A is stored; = ‘L’: Lower triangle of A is stored.

in
n

The number of linear equations, i.e., 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.

inout
A

Single complex array, dimension (lda, n). On entry, the Hermitian matrix A. If uplo = ‘U’, 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 = ‘L’, 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.

in
lda

The leading dimension of the array A. lda >= max(1, n).

out
ipiv

Integer array, dimension (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].

inout
B

Single complex array, dimension (ldb, nrhs). On entry, the N-by-NRHS right hand side matrix B. On exit, if info = 0, the N-by-NRHS solution matrix X.

in
ldb

The leading dimension of the array B. ldb >= max(1, n).

out
work

Single complex array, dimension (max(1, lwork)). On exit, if info = 0, work[0] returns the optimal lwork.

in
lwork

The length of work. lwork >= max(1, 2*n, 3*n-2). If lwork = -1, then a workspace query is assumed.

out
info

  • = 0: successful exit

  • < 0: if info = -i, the i-th argument had an illegal value

  • > 0: if info = i, D(i,i) is exactly zero.

Functions

void zhesv_aa(
    const char*          uplo,
    const INT            n,
    const INT            nrhs,
          c128* restrict A,
    const INT            lda,
          INT*  restrict ipiv,
          c128* restrict B,
    const INT            ldb,
          c128* restrict work,
    const INT            lwork,
          INT*           info
);
void zhesv_aa(const char *uplo, const INT n, const INT nrhs, c128 *restrict A, const INT lda, INT *restrict ipiv, c128 *restrict B, const INT ldb, c128 *restrict work, const INT lwork, INT *info)#

ZHESV_AA computes the solution to a complex system of linear equations A * X = B, where A is an N-by-N Hermitian matrix and X and B are N-by-NRHS matrices.

Aasen’s algorithm is used to factor A as A = U**H * T * U, if UPLO = ‘U’, or A = L * T * L**H, if UPLO = ‘L’, where U (or L) is a product of permutation and unit upper (lower) triangular matrices, and T is Hermitian and tridiagonal. The factored form of A is then used to solve the system of equations A * X = B.

On exit, if info = 0, the tridiagonal matrix T and the multipliers used to obtain the factor U or L from the factorization A = U**H*T*U or A = L*T*L**H as computed by ZHETRF_AA.

Parameters

in
uplo

= ‘U’: Upper triangle of A is stored; = ‘L’: Lower triangle of A is stored.

in
n

The number of linear equations, i.e., 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.

inout
A

Double complex array, dimension (lda, n). On entry, the Hermitian matrix A. If uplo = ‘U’, 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 = ‘L’, 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.

in
lda

The leading dimension of the array A. lda >= max(1, n).

out
ipiv

Integer array, dimension (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].

inout
B

Double complex array, dimension (ldb, nrhs). On entry, the N-by-NRHS right hand side matrix B. On exit, if info = 0, the N-by-NRHS solution matrix X.

in
ldb

The leading dimension of the array B. ldb >= max(1, n).

out
work

Double complex array, dimension (max(1, lwork)). On exit, if info = 0, work[0] returns the optimal lwork.

in
lwork

The length of work. lwork >= max(1, 2*n, 3*n-2). If lwork = -1, then a workspace query is assumed.

out
info

  • = 0: successful exit

  • < 0: if info = -i, the i-th argument had an illegal value

  • > 0: if info = i, D(i,i) is exactly zero.