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

Functions

int64_t lapack::pbcon (lapack::Uplo uplo, int64_t n, int64_t kd, double const *AB, int64_t ldab, double anorm, double *rcond)
 
int64_t lapack::pbcon (lapack::Uplo uplo, int64_t n, int64_t kd, float const *AB, int64_t ldab, float anorm, float *rcond)
 
int64_t lapack::pbcon (lapack::Uplo uplo, int64_t n, int64_t kd, std::complex< double > const *AB, int64_t ldab, double anorm, double *rcond)
 Estimates the reciprocal of the condition number (in the 1-norm) of a Hermitian positive definite band matrix using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::pbtrf.
 
int64_t lapack::pbcon (lapack::Uplo uplo, int64_t n, int64_t kd, std::complex< float > const *AB, int64_t ldab, float anorm, float *rcond)
 
int64_t lapack::pbequ (lapack::Uplo uplo, int64_t n, int64_t kd, double const *AB, int64_t ldab, double *S, double *scond, double *amax)
 
int64_t lapack::pbequ (lapack::Uplo uplo, int64_t n, int64_t kd, float const *AB, int64_t ldab, float *S, float *scond, float *amax)
 
int64_t lapack::pbequ (lapack::Uplo uplo, int64_t n, int64_t kd, std::complex< double > const *AB, int64_t ldab, double *S, double *scond, double *amax)
 Computes row and column scalings intended to equilibrate a Hermitian positive definite band matrix A and reduce its condition number (with respect to the two-norm).
 
int64_t lapack::pbequ (lapack::Uplo uplo, int64_t n, int64_t kd, std::complex< float > const *AB, int64_t ldab, float *S, float *scond, float *amax)
 
int64_t lapack::pbrfs (lapack::Uplo uplo, int64_t n, int64_t kd, int64_t nrhs, double const *AB, int64_t ldab, double const *AFB, int64_t ldafb, double const *B, int64_t ldb, double *X, int64_t ldx, double *ferr, double *berr)
 
int64_t lapack::pbrfs (lapack::Uplo uplo, int64_t n, int64_t kd, int64_t nrhs, float const *AB, int64_t ldab, float const *AFB, int64_t ldafb, float const *B, int64_t ldb, float *X, int64_t ldx, float *ferr, float *berr)
 
int64_t lapack::pbrfs (lapack::Uplo uplo, int64_t n, int64_t kd, int64_t nrhs, std::complex< double > const *AB, int64_t ldab, std::complex< double > const *AFB, int64_t ldafb, 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 banded, and provides error bounds and backward error estimates for the solution.
 
int64_t lapack::pbrfs (lapack::Uplo uplo, int64_t n, int64_t kd, int64_t nrhs, std::complex< float > const *AB, int64_t ldab, std::complex< float > const *AFB, int64_t ldafb, std::complex< float > const *B, int64_t ldb, std::complex< float > *X, int64_t ldx, float *ferr, float *berr)
 
int64_t lapack::pbstf (lapack::Uplo uplo, int64_t n, int64_t kd, double *AB, int64_t ldab)
 
int64_t lapack::pbstf (lapack::Uplo uplo, int64_t n, int64_t kd, float *AB, int64_t ldab)
 
int64_t lapack::pbstf (lapack::Uplo uplo, int64_t n, int64_t kd, std::complex< double > *AB, int64_t ldab)
 Computes a split Cholesky factorization of a Hermitian positive definite band matrix A.
 
int64_t lapack::pbstf (lapack::Uplo uplo, int64_t n, int64_t kd, std::complex< float > *AB, int64_t ldab)
 
int64_t lapack::pbtrf (lapack::Uplo uplo, int64_t n, int64_t kd, double *AB, int64_t ldab)
 
int64_t lapack::pbtrf (lapack::Uplo uplo, int64_t n, int64_t kd, float *AB, int64_t ldab)
 
int64_t lapack::pbtrf (lapack::Uplo uplo, int64_t n, int64_t kd, std::complex< double > *AB, int64_t ldab)
 Computes the Cholesky factorization of a Hermitian positive definite band matrix A.
 
int64_t lapack::pbtrf (lapack::Uplo uplo, int64_t n, int64_t kd, std::complex< float > *AB, int64_t ldab)
 
int64_t lapack::pbtrs (lapack::Uplo uplo, int64_t n, int64_t kd, int64_t nrhs, double const *AB, int64_t ldab, double *B, int64_t ldb)
 
int64_t lapack::pbtrs (lapack::Uplo uplo, int64_t n, int64_t kd, int64_t nrhs, float const *AB, int64_t ldab, float *B, int64_t ldb)
 
int64_t lapack::pbtrs (lapack::Uplo uplo, int64_t n, int64_t kd, int64_t nrhs, std::complex< double > const *AB, int64_t ldab, std::complex< double > *B, int64_t ldb)
 Solves a system of linear equations \(A X = B\) with a Hermitian positive definite band matrix A using the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) computed by lapack::pbtrf.
 
int64_t lapack::pbtrs (lapack::Uplo uplo, int64_t n, int64_t kd, int64_t nrhs, std::complex< float > const *AB, int64_t ldab, std::complex< float > *B, int64_t ldb)
 

Detailed Description

Function Documentation

◆ pbcon()

int64_t lapack::pbcon ( lapack::Uplo  uplo,
int64_t  n,
int64_t  kd,
std::complex< double > const *  AB,
int64_t  ldab,
double  anorm,
double *  rcond 
)

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

An estimate is obtained for \(|| A^{-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 triangular factor stored in AB;
  • lapack::Uplo::Lower: Lower triangular factor stored in AB.
[in]nThe order of the matrix A. n >= 0.
[in]kd
  • If uplo = Upper, the number of superdiagonals of the matrix A;
  • if uplo = Lower, the number of subdiagonals.
  • kd >= 0.
[in]ABThe n-by-n band matrix AB, stored in an ldab-by-n array. The triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) of the band matrix A, stored in the first kd+1 rows of the array. The j-th column of U or L is stored in the j-th column of the array AB as follows:
  • if uplo = Upper, AB(kd+1+i-j,j) = U(i,j) for max(1,j-kd) <= i <= j;
  • if uplo = Lower, AB(1+i-j,j) = L(i,j) for j <= i <= min(n,j+kd).
[in]ldabThe leading dimension of the array AB. ldab >= kd+1.
[in]anormThe 1-norm (or infinity-norm) of the Hermitian band 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 inv(A) computed in this routine.
Returns
= 0: successful exit

◆ pbequ()

int64_t lapack::pbequ ( lapack::Uplo  uplo,
int64_t  n,
int64_t  kd,
std::complex< double > const *  AB,
int64_t  ldab,
double *  S,
double *  scond,
double *  amax 
)

Computes row and column scalings intended to equilibrate a Hermitian positive definite band matrix A 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 triangular of A is stored;
  • lapack::Uplo::Lower: Lower triangular of A is stored.
[in]nThe order of the matrix A. n >= 0.
[in]kd
  • If uplo = Upper, the number of superdiagonals of the matrix A;
  • if uplo = Lower, the number of subdiagonals.
  • kd >= 0.
[in]ABThe n-by-n band matrix AB, stored in an ldab-by-n array. The upper or lower triangle of the Hermitian band matrix A, stored in the first kd+1 rows of the array. The j-th column of A is stored in the j-th column of the array AB as follows:
  • if uplo = Upper, AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd) <= i <= j;
  • if uplo = Lower, AB(1+i-j,j) = A(i,j) for j <= i <= min(n,j+kd).
[in]ldabThe leading dimension of the array A. ldab >= kd+1.
[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.

◆ pbrfs()

int64_t lapack::pbrfs ( lapack::Uplo  uplo,
int64_t  n,
int64_t  kd,
int64_t  nrhs,
std::complex< double > const *  AB,
int64_t  ldab,
std::complex< double > const *  AFB,
int64_t  ldafb,
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 banded, 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]kd
  • If uplo = Upper, the number of superdiagonals of the matrix A;
  • if uplo = Lower, the number of subdiagonals.
  • kd >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrices B and X. nrhs >= 0.
[in]ABThe n-by-n band matrix AB, stored in an ldab-by-n array. The upper or lower triangle of the Hermitian band matrix A, stored in the first kd+1 rows of the array. The j-th column of A is stored in the j-th column of the array AB as follows:
  • if uplo = Upper, AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd) <= i <= j;
  • if uplo = Lower, AB(1+i-j,j) = A(i,j) for j <= i <= min(n,j+kd).
[in]ldabThe leading dimension of the array AB. ldab >= kd+1.
[in]AFBThe n-by-n band matrix AFB, stored in an ldafb-by-n array. The triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) of the band matrix A as computed by lapack::pbtrf, in the same storage format as A (see AB).
[in]ldafbThe leading dimension of the array AFB. ldafb >= kd+1.
[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::pbtrs. 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

◆ pbstf()

int64_t lapack::pbstf ( lapack::Uplo  uplo,
int64_t  n,
int64_t  kd,
std::complex< double > *  AB,
int64_t  ldab 
)

Computes a split Cholesky factorization of a Hermitian positive definite band matrix A.

This routine is designed to be used in conjunction with lapack::hbgst.

The factorization has the form \(A = S^H S\) where S is a band matrix of the same bandwidth as A and the following structure:

\[ S = \begin{bmatrix} U & \\ M & L \end{bmatrix}, \]

where U is upper triangular of order m = (n+kd)/2, and L is lower triangular of order n-m.

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]kd
  • If uplo = Upper, the number of superdiagonals of the matrix A;
  • if uplo = Lower, the number of subdiagonals.
  • kd >= 0.
[in,out]ABThe n-by-n band matrix AB, stored in an ldab-by-n array.
  • On entry, the upper or lower triangle of the Hermitian band matrix A, stored in the first kd+1 rows of the array. The j-th column of A is stored in the j-th column of the array AB as follows:
    • if uplo = Upper, AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd) <= i <= j;
    • if uplo = Lower, AB(1+i-j,j) = A(i,j) for j <= i <= min(n,j+kd).
  • On successful exit, the factor S from the split Cholesky factorization \(A = S^H S.\) See Further Details.
[in]ldabThe leading dimension of the array AB. ldab >= kd+1.
Returns
= 0: successful exit
> 0: if return value = i, the factorization could not be completed, because the updated element a(i,i) was negative; the matrix A is not positive definite.
Further Details

The band storage scheme is illustrated by the following example, when n = 7, kd = 2:

S = [ s11  s12  s13                     ]
    [      s22  s23  s24                ]
    [           s33  s34                ]
    [                s44                ]
    [           s53  s54  s55           ]
    [                s64  s65  s66      ]
    [                     s75  s76  s77 ]

If uplo = Upper, the array AB holds:

on entry:

[  *    *   a13  a24  a35  a46  a57 ]
[  *   a12  a23  a34  a45  a56  a67 ]
[ a11  a22  a33  a44  a55  a66  a77 ]

on exit:

[  *    *   s13  s24  s53^H  s64^H  s75^H ]
[  *   s12  s23  s34  s54^H  s65^H  s76^H ]
[ s11  s22  s33  s44  s55    s66    s77   ]

If uplo = Lower, the array AB holds:

on entry:

a11  a22  a33  a44  a55  a66  a77
a21  a32  a43  a54  a65  a76   *
a31  a42  a53  a64  a64   *    *

on exit:

s11    s22    s33    s44  s55  s66  s77
s12^H  s23^H  s34^H  s54  s65  s76   *
s13^H  s24^H  s53    s64  s75   *    *

Array elements marked * are not used by the routine; s12^H denotes conj(s12); the diagonal elements of S are real.

◆ pbtrf()

int64_t lapack::pbtrf ( lapack::Uplo  uplo,
int64_t  n,
int64_t  kd,
std::complex< double > *  AB,
int64_t  ldab 
)

Computes the Cholesky factorization of a Hermitian positive definite band matrix A.

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]kd
  • If uplo = Upper, the number of superdiagonals of the matrix A;
  • if uplo = Lower, the number of subdiagonals.
  • kd >= 0.
[in,out]ABThe n-by-n band matrix AB, stored in an ldab-by-n array.
  • On entry, the upper or lower triangle of the Hermitian band matrix A, stored in the first kd+1 rows of the array. The j-th column of A is stored in the j-th column of the array AB as follows:
    • if uplo = Upper, AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd) <= i <= j;
    • if uplo = Lower, AB(1+i-j,j) = A(i,j) for j <= i <= min(n,j+kd).
  • On successful exit, the triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) of the band matrix A, in the same storage format as A.
[in]ldabThe leading dimension of the array AB. ldab >= kd+1.
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 band storage scheme is illustrated by the following example, when n = 6, kd = 2, and uplo = Upper:

On entry:                        On exit:

 *    *   a13  a24  a35  a46      *    *   u13  u24  u35  u46
 *   a12  a23  a34  a45  a56      *   u12  u23  u34  u45  u56
a11  a22  a33  a44  a55  a66     u11  u22  u33  u44  u55  u66

Similarly, if uplo = Lower the format of A is as follows:

On entry:                        On exit:

a11  a22  a33  a44  a55  a66     l11  l22  l33  l44  l55  l66
a21  a32  a43  a54  a65   *      l21  l32  l43  l54  l65   *
a31  a42  a53  a64   *    *      l31  l42  l53  l64   *    *

Array elements marked * are not used by the routine.

◆ pbtrs()

int64_t lapack::pbtrs ( lapack::Uplo  uplo,
int64_t  n,
int64_t  kd,
int64_t  nrhs,
std::complex< double > const *  AB,
int64_t  ldab,
std::complex< double > *  B,
int64_t  ldb 
)

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

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

Parameters
[in]uplo
  • lapack::Uplo::Upper: Upper triangular factor stored in AB;
  • lapack::Uplo::Lower: Lower triangular factor stored in AB.
[in]nThe order of the matrix A. n >= 0.
[in]kd
  • If uplo = Upper, the number of superdiagonals of the matrix A;
  • if uplo = Lower, the number of subdiagonals.
  • kd >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0.
[in]ABThe n-by-n band matrix AB, stored in an ldab-by-n array. The triangular factor U or L from the Cholesky factorization \(A = U^H U\) or \(A = L L^H\) of the band matrix A, stored in the first kd+1 rows of the array. The j-th column of U or L is stored in the j-th column of the array AB as follows:
  • if uplo = Upper, AB(kd+1+i-j,j) = U(i,j) for max(1,j-kd) <= i <= j;
  • if uplo = Lower, AB(1+i-j,j) = L(i,j) for j <= i <= min(n,j+kd).
[in]ldabThe leading dimension of the array AB. ldab >= kd+1.
[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