LAPACK++ 2024.05.31
LAPACK C++ API
Loading...
Searching...
No Matches
Positive definite: Cholesky: packed

Functions

int64_t lapack::ppcon (lapack::Uplo uplo, int64_t n, double const *AP, double anorm, double *rcond)
 
int64_t lapack::ppcon (lapack::Uplo uplo, int64_t n, float const *AP, float anorm, float *rcond)
 
int64_t lapack::ppcon (lapack::Uplo uplo, int64_t n, std::complex< double > const *AP, double anorm, double *rcond)
 Estimates the reciprocal of the condition number (in the 1-norm) of a Hermitian positive definite packed matrix using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::pptrf.
 
int64_t lapack::ppcon (lapack::Uplo uplo, int64_t n, std::complex< float > const *AP, float anorm, float *rcond)
 
int64_t lapack::ppequ (lapack::Uplo uplo, int64_t n, double const *AP, double *S, double *scond, double *amax)
 
int64_t lapack::ppequ (lapack::Uplo uplo, int64_t n, float const *AP, float *S, float *scond, float *amax)
 
int64_t lapack::ppequ (lapack::Uplo uplo, int64_t n, std::complex< double > const *AP, double *S, double *scond, double *amax)
 Computes row and column scalings intended to equilibrate a Hermitian positive definite matrix A in packed storage and reduce its condition number (with respect to the two-norm).
 
int64_t lapack::ppequ (lapack::Uplo uplo, int64_t n, std::complex< float > const *AP, float *S, float *scond, float *amax)
 
int64_t lapack::pprfs (lapack::Uplo uplo, int64_t n, int64_t nrhs, double const *AP, double const *AFP, double const *B, int64_t ldb, double *X, int64_t ldx, double *ferr, double *berr)
 
int64_t lapack::pprfs (lapack::Uplo uplo, int64_t n, int64_t nrhs, float const *AP, float const *AFP, float const *B, int64_t ldb, float *X, int64_t ldx, float *ferr, float *berr)
 
int64_t lapack::pprfs (lapack::Uplo uplo, int64_t n, int64_t nrhs, std::complex< double > const *AP, std::complex< double > const *AFP, std::complex< double > const *B, int64_t ldb, std::complex< double > *X, int64_t ldx, double *ferr, double *berr)
 Improves the computed solution to a system of linear equations when the coefficient matrix is Hermitian positive definite and packed, and provides error bounds and backward error estimates for the solution.
 
int64_t lapack::pprfs (lapack::Uplo uplo, int64_t n, int64_t nrhs, std::complex< float > const *AP, std::complex< float > const *AFP, std::complex< float > const *B, int64_t ldb, std::complex< float > *X, int64_t ldx, float *ferr, float *berr)
 
int64_t lapack::pptrf (lapack::Uplo uplo, int64_t n, double *AP)
 
int64_t lapack::pptrf (lapack::Uplo uplo, int64_t n, float *AP)
 
int64_t lapack::pptrf (lapack::Uplo uplo, int64_t n, std::complex< double > *AP)
 Computes the Cholesky factorization of a Hermitian positive definite matrix A stored in packed format.
 
int64_t lapack::pptrf (lapack::Uplo uplo, int64_t n, std::complex< float > *AP)
 
int64_t lapack::pptri (lapack::Uplo uplo, int64_t n, double *AP)
 
int64_t lapack::pptri (lapack::Uplo uplo, int64_t n, float *AP)
 
int64_t lapack::pptri (lapack::Uplo uplo, int64_t n, std::complex< double > *AP)
 Computes the inverse of a Hermitian positive definite matrix A using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::pptrf.
 
int64_t lapack::pptri (lapack::Uplo uplo, int64_t n, std::complex< float > *AP)
 
int64_t lapack::pptrs (lapack::Uplo uplo, int64_t n, int64_t nrhs, double const *AP, double *B, int64_t ldb)
 
int64_t lapack::pptrs (lapack::Uplo uplo, int64_t n, int64_t nrhs, float const *AP, float *B, int64_t ldb)
 
int64_t lapack::pptrs (lapack::Uplo uplo, int64_t n, int64_t nrhs, std::complex< double > const *AP, std::complex< double > *B, int64_t ldb)
 Solves a system of linear equations \(A X = B\) with a Hermitian positive definite matrix A in packed storage using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::pptrf.
 
int64_t lapack::pptrs (lapack::Uplo uplo, int64_t n, int64_t nrhs, std::complex< float > const *AP, std::complex< float > *B, int64_t ldb)
 

Detailed Description

Function Documentation

◆ ppcon()

int64_t lapack::ppcon ( lapack::Uplo  uplo,
int64_t  n,
std::complex< double > const *  AP,
double  anorm,
double *  rcond 
)

Estimates the reciprocal of the condition number (in the 1-norm) of a Hermitian positive definite packed matrix using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::pptrf.

An estimate is obtained for \(||A^{-1}||_1\), and the reciprocal of the condition number is computed as \(\text{rcond} = 1 / (||A||_1 \cdot ||A^{-1}||_1).\)

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 order of the matrix A. n >= 0.
[in]APThe vector AP of length n*(n+1)/2. The triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H,\) packed columnwise in a linear array. The j-th column of U or L is stored in the array AP as follows:
  • if uplo = Upper, AP(i + (j-1)*j/2) = U(i,j) for 1 <= i <= j;
  • if uplo = Lower, AP(i + (j-1)*(2n-j)/2) = L(i,j) for j <= i <= n.
[in]anormThe 1-norm (or infinity-norm) of the Hermitian 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

◆ ppequ()

int64_t lapack::ppequ ( lapack::Uplo  uplo,
int64_t  n,
std::complex< double > const *  AP,
double *  S,
double *  scond,
double *  amax 
)

Computes row and column scalings intended to equilibrate a Hermitian positive definite matrix A in packed storage and reduce its condition number (with respect to the two-norm).

S contains the scale factors, \(S_i = 1 / \sqrt{ A_{i,i} },\) chosen so that the scaled matrix B with elements \(B_{i,j} = S_{i} A_{i,j} S_{j}\) has ones on the diagonal. This choice of S puts the condition number of B within a factor n of the smallest possible condition number over all possible diagonal scalings.

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 order of the matrix A. n >= 0.
[in]APThe vector AP of length n*(n+1)/2. 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.
[out]SThe vector S of length n. If successful, S contains the scale factors for A.
[out]scondIf successful, S contains the ratio of the smallest S(i) to the largest S(i). If scond >= 0.1 and amax is neither too large nor too small, it is not worth scaling by S.
[out]amaxAbsolute value of largest matrix element. If amax is very close to overflow or very close to underflow, the matrix should be scaled.
Returns
= 0: successful exit
> 0: if return value = i, the i-th diagonal element is nonpositive.

◆ pprfs()

int64_t lapack::pprfs ( lapack::Uplo  uplo,
int64_t  n,
int64_t  nrhs,
std::complex< double > const *  AP,
std::complex< double > const *  AFP,
std::complex< double > const *  B,
int64_t  ldb,
std::complex< double > *  X,
int64_t  ldx,
double *  ferr,
double *  berr 
)

Improves the computed solution to a system of linear equations when the coefficient matrix is Hermitian positive definite and packed, and provides error bounds and backward error estimates for the solution.

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 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]APThe vector AP of length n*(n+1)/2. 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.
[in]AFPThe vector AFP of length n*(n+1)/2. The triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H,\) as computed by lapack::pptrf, packed columnwise in a linear array in the same format as A (see AP).
[in]BThe n-by-nrhs matrix B, stored in an ldb-by-nrhs array. The right hand side matrix B.
[in]ldbThe leading dimension of the array B. ldb >= max(1,n).
[in,out]XThe n-by-nrhs matrix X, stored in an ldx-by-nrhs array. On entry, the solution matrix X, as computed by lapack::pptrs. On exit, the improved solution matrix X.
[in]ldxThe leading dimension of the array X. ldx >= max(1,n).
[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

◆ pptrf()

int64_t lapack::pptrf ( lapack::Uplo  uplo,
int64_t  n,
std::complex< double > *  AP 
)

Computes the Cholesky factorization of a Hermitian positive definite matrix A stored in packed format.

The factorization has the form \(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 lower triangular.

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 order of the matrix A. n >= 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 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.
Returns
= 0: successful exit
> 0: if return value = i, the leading minor of order i is not positive definite, and the factorization could not be completed.
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 ]

◆ pptri()

int64_t lapack::pptri ( lapack::Uplo  uplo,
int64_t  n,
std::complex< double > *  AP 
)

Computes the inverse of a Hermitian positive definite matrix A using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::pptrf.

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

Parameters
[in]uplo
  • lapack::Uplo::Upper: Upper triangular factor is stored in AP;
  • lapack::Uplo::Lower: Lower triangular factor is stored in AP.
[in]nThe order of the matrix A. n >= 0.
[in,out]APThe vector AP of length n*(n+1)/2.
  • On entry, the triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H,\) packed columnwise as a linear array. The j-th column of U or L is stored in the array AP as follows:
    • if uplo = Upper, AP(i + (j-1)*j/2) = U(i,j) for 1 <= i <= j;
    • if uplo = Lower, AP(i + (j-1)*(2n-j)/2) = L(i,j) for j <= i <= n.
  • On exit, the upper or lower triangle of the (Hermitian) inverse of A, overwriting the input factor U or L.
Returns
= 0: successful exit
> 0: if return value = i, the (i,i) element of the factor U or L is zero, and the inverse could not be computed.

◆ pptrs()

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

Solves a system of linear equations \(A X = B\) with a Hermitian positive definite matrix A in packed storage using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::pptrf.

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 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]APThe vector AP of length n*(n+1)/2. The triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H,\) packed columnwise in a linear array. The j-th column of U or L is stored in the array AP as follows:
  • if uplo = Upper, AP(i + (j-1)*j/2) = U(i,j) for 1 <= i <= j;
  • if uplo = Lower, AP(i + (j-1)*(2n-j)/2) = L(i,j) for j <= i <= n.
[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