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

Functions

int64_t lapack::ptcon (int64_t n, double const *D, double const *E, double anorm, double *rcond)
 
int64_t lapack::ptcon (int64_t n, double const *D, std::complex< double > const *E, double anorm, double *rcond)
 Computes the reciprocal of the condition number (in the 1-norm) of a Hermitian positive definite tridiagonal matrix using the factorization \(A = L D L^H\) or \(A = U^H D U\) computed by lapack::pttrf.
 
int64_t lapack::ptcon (int64_t n, float const *D, float const *E, float anorm, float *rcond)
 
int64_t lapack::ptcon (int64_t n, float const *D, std::complex< float > const *E, float anorm, float *rcond)
 
int64_t lapack::ptrfs (int64_t n, int64_t nrhs, double const *D, double const *E, double const *DF, double const *EF, double const *B, int64_t ldb, double *X, int64_t ldx, double *ferr, double *berr)
 
int64_t lapack::ptrfs (int64_t n, int64_t nrhs, float const *D, float const *E, float const *DF, float const *EF, float const *B, int64_t ldb, float *X, int64_t ldx, float *ferr, float *berr)
 
int64_t lapack::ptrfs (lapack::Uplo uplo, int64_t n, int64_t nrhs, double const *D, std::complex< double > const *E, double const *DF, std::complex< double > const *EF, 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 tridiagonal, and provides error bounds and backward error estimates for the solution.
 
int64_t lapack::ptrfs (lapack::Uplo uplo, int64_t n, int64_t nrhs, float const *D, std::complex< float > const *E, float const *DF, std::complex< float > const *EF, std::complex< float > const *B, int64_t ldb, std::complex< float > *X, int64_t ldx, float *ferr, float *berr)
 
int64_t lapack::pttrf (int64_t n, double *D, double *E)
 
int64_t lapack::pttrf (int64_t n, double *D, std::complex< double > *E)
 Computes the \(L D L^H\) factorization of a Hermitian positive definite tridiagonal matrix A.
 
int64_t lapack::pttrf (int64_t n, float *D, float *E)
 
int64_t lapack::pttrf (int64_t n, float *D, std::complex< float > *E)
 
int64_t lapack::pttrs (int64_t n, int64_t nrhs, double const *D, double const *E, double *B, int64_t ldb)
 
int64_t lapack::pttrs (int64_t n, int64_t nrhs, float const *D, float const *E, float *B, int64_t ldb)
 
int64_t lapack::pttrs (lapack::Uplo uplo, int64_t n, int64_t nrhs, double const *D, std::complex< double > const *E, std::complex< double > *B, int64_t ldb)
 Solves a tridiagonal system of the form.
 
int64_t lapack::pttrs (lapack::Uplo uplo, int64_t n, int64_t nrhs, float const *D, std::complex< float > const *E, std::complex< float > *B, int64_t ldb)
 

Detailed Description

Function Documentation

◆ ptcon()

int64_t lapack::ptcon ( int64_t  n,
double const *  D,
std::complex< double > const *  E,
double  anorm,
double *  rcond 
)

Computes the reciprocal of the condition number (in the 1-norm) of a Hermitian positive definite tridiagonal matrix using the factorization \(A = L D L^H\) or \(A = U^H D U\) computed by lapack::pttrf.

\(|| A^{-1} ||\) is computed by a direct method, 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]nThe order of the matrix A. n >= 0.
[in]DThe vector D of length n. The n diagonal elements of the diagonal matrix D from the factorization of A, as computed by lapack::pttrf.
[in]EThe vector E of length n-1. The (n-1) off-diagonal elements of the unit bidiagonal factor U or L from the factorization of A, as computed by lapack::pttrf.
[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 the 1-norm of \(A^{-1}\) computed in this routine.
Returns
= 0: successful exit
Further Details

The method used is described in Nicholas J. Higham, "Efficient Algorithms for Computing the Condition Number of a Tridiagonal Matrix", SIAM J. Sci. Stat. Comput., Vol. 7, No. 1, January 1986.

◆ ptrfs()

int64_t lapack::ptrfs ( lapack::Uplo  uplo,
int64_t  n,
int64_t  nrhs,
double const *  D,
std::complex< double > const *  E,
double const *  DF,
std::complex< double > const *  EF,
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 tridiagonal, 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]uploWhether the superdiagonal or the subdiagonal of the tridiagonal matrix A is stored and the form of the factorization:
  • lapack::Uplo::Upper: E is the superdiagonal of A, and \(A = U^H D U;\)
  • lapack::Uplo::Lower: E is the subdiagonal of A, and \(A = L D L^H.\)
    (The two forms are equivalent if A is real.)
[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]DThe vector D of length n. The n real diagonal elements of the tridiagonal matrix A.
[in]EThe vector E of length n-1. The (n-1) off-diagonal elements of the tridiagonal matrix A (see uplo).
[in]DFThe vector DF of length n. The n diagonal elements of the diagonal matrix D from the factorization computed by lapack::pttrf.
[in]EFThe vector EF of length n-1. The (n-1) off-diagonal elements of the unit bidiagonal factor U or L from the factorization computed by lapack::pttrf (see uplo).
[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::pttrs. 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 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).
[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

◆ pttrf()

int64_t lapack::pttrf ( int64_t  n,
double *  D,
std::complex< double > *  E 
)

Computes the \(L D L^H\) factorization of a Hermitian positive definite tridiagonal matrix A.

The factorization may also be regarded as having the form \(A = U^H D U.\)

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

Parameters
[in]nThe order of the matrix A. n >= 0.
[in,out]DThe vector D of length n. On entry, the n diagonal elements of the tridiagonal matrix A. On exit, the n diagonal elements of the diagonal matrix D from the \(L D L^H\) factorization of A.
[in,out]EThe vector E of length n-1. On entry, the (n-1) subdiagonal elements of the tridiagonal matrix A. On exit, the (n-1) subdiagonal elements of the unit bidiagonal factor L from the \(L D L^H\) factorization of A. E can also be regarded as the superdiagonal of the unit bidiagonal factor U from the \(U^H D U\) factorization of A.
Returns
= 0: successful exit
> 0: if return value = i, the leading minor of order i is not positive definite; if i < n, the factorization could not be completed, while if i = n, the factorization was completed, but D(n) <= 0.

◆ pttrs()

int64_t lapack::pttrs ( lapack::Uplo  uplo,
int64_t  n,
int64_t  nrhs,
double const *  D,
std::complex< double > const *  E,
std::complex< double > *  B,
int64_t  ldb 
)

Solves a tridiagonal system of the form.

\[ A X = B \]

using the factorization \(A = U^H D U\) or \(A = L D L^H\) computed by lapack::pttrf. D is a diagonal matrix specified in the vector D, U (or L) is a unit bidiagonal matrix whose superdiagonal (subdiagonal) is specified in the vector E, and X and B are n by nrhs matrices.

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

Parameters
[in]uploSpecifies the form of the factorization and whether the vector E is the superdiagonal of the upper bidiagonal factor U or the subdiagonal of the lower bidiagonal factor L.
  • lapack::Uplo::Upper: \(A = U^H D U,\) E is the superdiagonal of U
  • lapack::Uplo::Lower: \(A = L D L^H,\) E is the subdiagonal of L
[in]nThe order of the tridiagonal matrix A. n >= 0.
[in]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 0.
[in]DThe vector D of length n. The n diagonal elements of the diagonal matrix D from the factorization \(A = U^H D U\) or \(A = L D L^H.\)
[in]EThe vector E of length n-1.
  • If uplo = Upper, the (n-1) superdiagonal elements of the unit bidiagonal factor U from the factorization \(A = U^H D U.\)
  • If uplo = Lower, the (n-1) subdiagonal elements of the unit bidiagonal factor L from the factorization \(A = L D L^H.\)
[in,out]BThe n-by-nrhs matrix B, stored in an ldb-by-nrhs array. On entry, the right hand side vectors B for the system of linear equations. On exit, the solution vectors, X.
[in]ldbThe leading dimension of the array B. ldb >= max(1,n).
Returns
= 0: successful exit