LAPACK++ 2024.05.31
LAPACK C++ API
Loading...
Searching...
No Matches
Symmetric indefinite: Rook

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)
 

Detailed Description

Macro Definition Documentation

◆ sytrf_3

#define sytrf_3   sytrf_rk
See also
lapack::sytrf_rk

◆ sytri_3

#define sytri_3   sytri_rk
See also
lapack::sytri_rk

◆ sytrs_3

#define sytrs_3   sytrs_rk
See also
lapack::sytrs_rk

Function Documentation

◆ sycon_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.

Since
LAPACK 3.7.0.
Parameters
[in]uploSpecifies whether the details of the factorization are stored as an upper or lower triangular matrix:
  • lapack::Uplo::Upper: Upper triangular, form is \(A = P U D U^T P^T;\)
  • lapack::Uplo::Lower: Lower triangular, form is \(A = P L D L^T P^T.\)
[in]nThe order of the matrix A. n >= 0.
[in]AThe 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:
  • ONLY diagonal elements of the symmetric block diagonal matrix D on the diagonal of A, i.e. D(k,k) = A(k,k); (superdiagonal (or subdiagonal) elements of D should be provided on entry in array E), and
  • If uplo = Upper: factor U in the superdiagonal part of A.
  • If uplo = Lower: factor L in the subdiagonal part of A.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in]EThe 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
  • If uplo = Upper: E(i) = D(i-1,i),i=2:n, E(1) not referenced;
  • If uplo = Lower: E(i) = D(i+1,i),i=1:n-1, E(n) not referenced.
  • Note: For 1-by-1 diagonal block D(k), where 1 <= k <= n, the element E(k) is not referenced in both uplo = Upper or uplo = Lower cases.
[in]ipivThe vector ipiv of length n. Details of the interchanges and the block structure of D as determined by lapack::sytrf_rk.
[in]anormThe 1-norm of the original matrix A.
[out]rcondThe 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.
Returns
= 0: successful exit

◆ sytrf_rk()

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.

Since
LAPACK 3.7.0. This interface replaces the older lapack::sytrf_rook.
Parameters
[in]uploWhether the upper or lower triangular part of the symmetric matrix A is stored:
  • lapack::Uplo::Upper: Upper triangular
  • lapack::Uplo::Lower: Lower triangular
[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 symmetric 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, contains:
    • ONLY diagonal elements of the symmetric block diagonal matrix D on the diagonal of A, i.e. D(k,k) = A(k,k); (superdiagonal (or subdiagonal) elements of D are stored on exit in array E), and
    • If uplo = Upper: factor U in the superdiagonal part of A.
    • If uplo = Lower: factor L in the subdiagonal part of A.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[out]EThe 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
  • If uplo = Upper: E(i) = D(i-1,i), i=2:n, E(1) is set to 0;
  • If uplo = Lower: E(i) = D(i+1,i), i=1:n-1, E(n) is set to 0.
  • Note: For 1-by-1 diagonal block D(k), where 1 <= k <= n, the element E(k) is set to 0 in both uplo = Upper or uplo = Lower cases.
[out]ipivThe 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.
  • If uplo = Upper, (in factorization order, k decreases from n to 1): a) A single positive entry ipiv(k) > 0 means: D(k,k) is a 1-by-1 diagonal block. If ipiv(k) != k, rows and columns k and ipiv(k) were interchanged in the matrix A(1:n,1:n); If ipiv(k) = k, no interchange occurred.

    b) A pair of consecutive negative entries ipiv(k) < 0 and ipiv(k-1) < 0 means: D(k-1:k,k-1:k) is a 2-by-2 diagonal block. (NOTE: negative entries in ipiv appear ONLY in pairs). 1) If -ipiv(k) != k, rows and columns k and -ipiv(k) were interchanged in the matrix A(1:n,1:n). If -ipiv(k) = k, no interchange occurred. 2) If -ipiv(k-1) != k-1, rows and columns k-1 and -ipiv(k-1) were interchanged in the matrix A(1:n,1:n). If -ipiv(k-1) = k-1, no interchange occurred.

    c) In both cases a) and b), always ABS( ipiv(k) ) <= k.

    d) NOTE: Any entry ipiv(k) is always NONZERO on output.

  • If uplo = Lower, (in factorization order, k increases from 1 to n): a) A single positive entry ipiv(k) > 0 means: D(k,k) is a 1-by-1 diagonal block. If ipiv(k) != k, rows and columns k and ipiv(k) were interchanged in the matrix A(1:n,1:n). If ipiv(k) = k, no interchange occurred.

    b) A pair of consecutive negative entries ipiv(k) < 0 and ipiv(k+1) < 0 means: D(k:k+1,k:k+1) is a 2-by-2 diagonal block. (NOTE: negative entries in ipiv appear ONLY in pairs). 1) If -ipiv(k) != k, rows and columns k and -ipiv(k) were interchanged in the matrix A(1:n,1:n). If -ipiv(k) = k, no interchange occurred. 2) If -ipiv(k+1) != k+1, rows and columns k-1 and -ipiv(k-1) were interchanged in the matrix A(1:n,1:n). If -ipiv(k+1) = k+1, no interchange occurred.

    c) In both cases a) and b), always ABS( ipiv(k) ) >= k.

    d) NOTE: Any entry ipiv(k) is always NONZERO on output.

Returns
= 0: successful exit
> 0: If return value = i, the matrix A is singular, because: If uplo = Upper: column i in the upper triangular part of A contains all zeros. If uplo = Lower: column i in the lower triangular part of A contains all zeros.

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.

◆ sytrf_rook()

int64_t lapack::sytrf_rook ( lapack::Uplo  uplo,
int64_t  n,
std::complex< double > *  A,
int64_t  lda,
int64_t *  ipiv 
)
See also
lapack::sytrf_rk.
Since
LAPACK 3.5.0.

◆ sytri_rk()

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.

Since
LAPACK 3.7.0. This interface replaces the older lapack::sytri_rook.
Parameters
[in]uploSpecifies whether the details of the factorization are stored as an upper or lower triangular matrix.
  • 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, diagonal of the block diagonal matrix D and factors U or L as computed by lapack::sytrf_rk:
    • ONLY diagonal elements of the symmetric block diagonal matrix D on the diagonal of A, i.e. D(k,k) = A(k,k); (superdiagonal (or subdiagonal) elements of D should be provided on entry in array E), and
    • If uplo = Upper: factor U in the superdiagonal part of A.
    • If uplo = Lower: factor L in the subdiagonal part of A.
  • On successful exit, the symmetric inverse of the original matrix.
    • If uplo = Upper: the upper triangular part of the inverse is formed and the part of A below the diagonal is not referenced;
    • If uplo = Lower: the lower triangular part of the inverse is formed and the part of A above the diagonal is not referenced.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in]EThe 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
  • If uplo = Upper: E(i) = D(i-1,i),i=2:n, E(1) not referenced;
  • If uplo = Lower: E(i) = D(i+1,i),i=1:n-1, E(n) not referenced.
  • Note: For 1-by-1 diagonal block D(k), where 1 <= k <= n, the element E(k) is not referenced in both uplo = Upper or uplo = Lower cases.
[in]ipivThe vector ipiv of length n. Details of the interchanges and the block structure of D as determined by lapack::sytrf_rk.
Returns
= 0: successful exit
> 0: if return value = i, D(i,i) = 0; the matrix is singular and its inverse could not be computed.

◆ sytrs_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.

Since
LAPACK 3.7.0. This interface replaces the older lapack::sytrs_rook.
Parameters
[in]uploSpecifies whether the details of the factorization are stored as an upper or lower triangular matrix:
  • lapack::Uplo::Upper: Upper triangular, form is \(A = P U D (U^T) (P^T);\)
  • lapack::Uplo::Lower: Lower triangular, form is \(A = P L D (L^T) (P^T).\)
[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. Diagonal of the block diagonal matrix D and factors U or L as computed by lapack::sytrf_rk:
  • ONLY diagonal elements of the symmetric block diagonal matrix D on the diagonal of A, i.e. D(k,k) = A(k,k); (superdiagonal (or subdiagonal) elements of D should be provided on entry in array E), and
  • If uplo = Upper: factor U in the superdiagonal part of A.
  • If uplo = Lower: factor L in the subdiagonal part of A.
[in]ldaThe leading dimension of the array A. lda >= max(1,n).
[in]EThe 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
  • If uplo = Upper: E(i) = D(i-1,i),i=2:n, E(1) not referenced;
  • If uplo = Lower: E(i) = D(i+1,i),i=1:n-1, E(n) not referenced.
  • Note: For 1-by-1 diagonal block D(k), where 1 <= k <= n, the element E(k) is not referenced in both uplo = Upper or uplo = Lower cases.
[in]ipivThe vector ipiv of length n. Details of the interchanges and the block structure of D as determined by lapack::sytrf_rk.
[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

◆ sytrs_rook()

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 
)
See also
lapack::sytrs_rk
Since
LAPACK 3.5.0.