LAPACK++  2022.07.00
LAPACK C++ API
Positive definite: Cholesky: packed

Functions

int64_t lapack::ppsv (lapack::Uplo uplo, int64_t n, int64_t nrhs, double *AP, double *B, int64_t ldb)
 
int64_t lapack::ppsv (lapack::Uplo uplo, int64_t n, int64_t nrhs, float *AP, float *B, int64_t ldb)
 
int64_t lapack::ppsv (lapack::Uplo uplo, int64_t n, int64_t nrhs, std::complex< double > *AP, std::complex< double > *B, int64_t ldb)
 Computes the solution to a system of linear equations. More...
 
int64_t lapack::ppsv (lapack::Uplo uplo, int64_t n, int64_t nrhs, std::complex< float > *AP, std::complex< float > *B, int64_t ldb)
 
int64_t lapack::ppsvx (lapack::Factored fact, lapack::Uplo uplo, int64_t n, int64_t nrhs, double *AP, double *AFP, lapack::Equed *equed, double *S, double *B, int64_t ldb, double *X, int64_t ldx, double *rcond, double *ferr, double *berr)
 
int64_t lapack::ppsvx (lapack::Factored fact, lapack::Uplo uplo, int64_t n, int64_t nrhs, float *AP, float *AFP, lapack::Equed *equed, float *S, float *B, int64_t ldb, float *X, int64_t ldx, float *rcond, float *ferr, float *berr)
 
int64_t lapack::ppsvx (lapack::Factored fact, lapack::Uplo uplo, int64_t n, int64_t nrhs, std::complex< double > *AP, std::complex< double > *AFP, lapack::Equed *equed, double *S, std::complex< double > *B, int64_t ldb, std::complex< double > *X, int64_t ldx, double *rcond, double *ferr, double *berr)
 Uses the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) to compute the solution to a system of linear equations. More...
 
int64_t lapack::ppsvx (lapack::Factored fact, lapack::Uplo uplo, int64_t n, int64_t nrhs, std::complex< float > *AP, std::complex< float > *AFP, lapack::Equed *equed, float *S, std::complex< float > *B, int64_t ldb, std::complex< float > *X, int64_t ldx, float *rcond, float *ferr, float *berr)
 

Detailed Description

Function Documentation

◆ ppsv()

int64_t lapack::ppsv ( lapack::Uplo  uplo,
int64_t  n,
int64_t  nrhs,
std::complex< double > *  AP,
std::complex< double > *  B,
int64_t  ldb 
)

Computes the solution to a system of linear equations.

\[ A X = B, \]

where A is an n-by-n Hermitian positive definite matrix stored in packed format and X and B are n-by-nrhs matrices.

The Cholesky decomposition is used to factor A as \(A = U^H U\) if uplo = Upper, or \(A = L L^H\) if uplo = Lower, where U is an upper triangular matrix and L is a lower triangular matrix. The factored form of A is then used to solve the system of equations \(A X = B.\)

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Parameters
[in]uplo
  • lapack::Uplo::Upper: Upper triangle of A is stored;
  • lapack::Uplo::Lower: Lower triangle of A is stored.
[in]nThe number of linear equations, i.e., the 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,out]APThe vector AP of length n*(n+1)/2.
  • On entry, the upper or lower triangle of the Hermitian matrix A, packed columnwise in a linear array. The j-th column of A is stored in the array AP as follows:
    • if uplo = Upper, AP(i + (j-1)*j/2) = A(i,j) for 1 <= i <= j;
    • if uplo = Lower, AP(i + (j-1)*(2n-j)/2) = A(i,j) for j <= i <= n.
      See below for further details.
  • On successful exit, the factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H,\) in the same storage format as A.
[in,out]BThe n-by-nrhs matrix B, stored in an ldb-by-nrhs array. On entry, the n-by-nrhs right hand side matrix B. On successful exit, the n-by-nrhs solution matrix X.
[in]ldbThe leading dimension of the array B. ldb >= max(1,n).
Returns
= 0: successful exit
> 0: if return value = i, the leading minor of order i of A is not positive definite, so the factorization could not be completed, and the solution has not been computed.
Further Details

The packed storage scheme is illustrated by the following example when n = 4, uplo = Upper:

Two-dimensional storage of the Hermitian matrix A:

[ a11 a12 a13 a14 ]
[     a22 a23 a24 ]
[         a33 a34 ]    (aij = conj(aji))
[             a44 ]

Packed storage of the upper triangle of A:

AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]

◆ ppsvx()

int64_t lapack::ppsvx ( lapack::Factored  fact,
lapack::Uplo  uplo,
int64_t  n,
int64_t  nrhs,
std::complex< double > *  AP,
std::complex< double > *  AFP,
lapack::Equed *  equed,
double *  S,
std::complex< double > *  B,
int64_t  ldb,
std::complex< double > *  X,
int64_t  ldx,
double *  rcond,
double *  ferr,
double *  berr 
)

Uses the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) to compute the solution to a system of linear equations.

\[ A X = B, \]

where A is an n-by-n Hermitian positive definite matrix stored in packed format and X and B are n-by-nrhs matrices.

Error bounds on the solution and a condition estimate are also provided.

Overloaded versions are available for float, double, std::complex<float>, and std::complex<double>.

Parameters
[in]factWhether or not the factored form of the matrix A is supplied on entry, and if not, whether the matrix A should be equilibrated before it is factored.
  • lapack::Factored::Factored: On entry, AFP contains the factored form of A. If equed = Yes, the matrix A has been equilibrated with scaling factors given by S. AP and AFP will not be modified.
  • lapack::Factored::NotFactored: The matrix A will be copied to AFP and factored.
  • lapack::Factored::Equilibrate: The matrix A will be equilibrated if necessary, then copied to AFP and factored.
[in]uplo
  • lapack::Uplo::Upper: Upper triangle of A is stored;
  • lapack::Uplo::Lower: Lower triangle of A is stored.
[in]nThe number of linear equations, i.e., the order of the matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrices B and X. nrhs >= 0.
[in,out]APThe vector AP of length n*(n+1)/2.
  • On entry, the upper or lower triangle of the Hermitian matrix A, packed columnwise in a linear array, except if fact = Factored and equed = Yes, then A must contain the equilibrated matrix \(\text{diag}(S) \; A \; \text{diag}(S)\). The j-th column of A is stored in the array AP as follows:
    • if uplo = Upper, AP(i + (j-1)*j/2) = A(i,j) for 1 <= i <= j;
    • if uplo = Lower, AP(i + (j-1)*(2n-j)/2) = A(i,j) for j <= i <= n.
      See below for further details.
  • A is not modified if fact = Factored or NotFactored, or if fact = Equilibrate and equed = None on exit.
  • On exit, if fact = Equilibrate and equed = Yes, A is overwritten by \(\text{diag}(S) \; A \; \text{diag}(S)\).
[in,out]AFPThe vector AFP of length n*(n+1)/2.
  • If fact = Factored, then AFP is an input argument and on entry contains the triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H,\) in the same storage format as A.
  • If equed != None, then AFP is the factored form of the equilibrated matrix A.
  • If fact = NotFactored, then AFP is an output argument and on exit returns the triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) of the original matrix A.
  • If fact = Equilibrate, then AFP is an output argument and on exit returns the triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) of the equilibrated matrix A (see the description of AP for the form of the equilibrated matrix).
[in,out]equedThe form of equilibration that was done:
  • lapack::Equed::None: No equilibration (always true if fact = NotFactored).
  • lapack::Equed::Yes: Equilibration was done, i.e., A has been replaced by \(\text{diag}(S) \; A \; \text{diag}(S).\)
  • equed is an input argument if fact = Factored; otherwise, it is an output argument.
[in,out]SThe vector S of length n. The scale factors for A.
  • If fact = Factored, S is an input argument;
  • otherwise, S is an output argument.
  • If fact = Factored and equed = Yes, each element of S must be positive.
  • If equed = None, S is not accessed.
[in,out]BThe n-by-nrhs matrix B, stored in an ldb-by-nrhs array. On entry, the n-by-nrhs right hand side matrix B. On exit, if equed = None, B is not modified; if equed = Yes, B is overwritten by \(\text{diag}(S) \; B.\)
[in]ldbThe leading dimension of the array B. ldb >= max(1,n).
[out]XThe n-by-nrhs matrix X, stored in an ldx-by-nrhs array. If successful or return value = n+1, the n-by-nrhs solution matrix X to the original system of equations. Note that if equed = Yes, A and B are modified on exit, and the solution to the equilibrated system is \(\text{diag}(S)^{-1} X.\)
[in]ldxThe leading dimension of the array X. ldx >= max(1,n).
[out]rcondThe estimate of the reciprocal condition number of the matrix A after equilibration (if done). If rcond is less than the machine precision (in particular, if rcond = 0), the matrix is singular to working precision. This condition is indicated by a return code of return value > 0.
[out]ferrThe vector ferr of length nrhs. The estimated forward error bound for each solution vector X(j) (the j-th column of the solution matrix X). If XTRUE is the true solution corresponding to X(j), ferr(j) is an estimated upper bound for the magnitude of the largest element in (X(j) - XTRUE) divided by the magnitude of the largest element in X(j). The estimate is as reliable as the estimate for rcond, and is almost always a slight overestimate of the true error.
[out]berrThe vector berr of length nrhs. The componentwise relative backward error of each solution vector X(j) (i.e., the smallest relative change in any element of A or B that makes X(j) an exact solution).
Returns
= 0: successful exit
> 0 and <= n: if return value = i, the leading minor of order i of A is not positive definite, so the factorization could not be completed, and the solution has not been computed. rcond = 0 is returned.
= n+1: U is nonsingular, but rcond is less than machine precision, meaning that the matrix is singular to working precision. Nevertheless, the solution and error bounds are computed because there are a number of situations where the computed solution can be more accurate than the value of rcond would suggest.
Further Details

The packed storage scheme is illustrated by the following example when n = 4, uplo = Upper:

Two-dimensional storage of the Hermitian matrix A:

[ a11 a12 a13 a14 ]
[     a22 a23 a24 ]
[         a33 a34 ]    (aij = conj(aji))
[             a44 ]

Packed storage of the upper triangle of A:

AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]