LAPACK++
LAPACK C++ API
Positive definite: Cholesky: tridiagonal

Functions

int64_t lapack::ptsv (int64_t n, int64_t nrhs, double *D, double *E, double *B, int64_t ldb)
 
int64_t lapack::ptsv (int64_t n, int64_t nrhs, double *D, std::complex< double > *E, 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 tridiagonal matrix, and X and B are n-by-nrhs matrices. More...
 
int64_t lapack::ptsv (int64_t n, int64_t nrhs, float *D, float *E, float *B, int64_t ldb)
 
int64_t lapack::ptsv (int64_t n, int64_t nrhs, float *D, std::complex< float > *E, std::complex< float > *B, int64_t ldb)
 
int64_t lapack::ptsvx (lapack::Factored fact, int64_t n, int64_t nrhs, double const *D, std::complex< double > const *E, double *DF, std::complex< double > *EF, std::complex< double > const *B, int64_t ldb, std::complex< double > *X, int64_t ldx, double *rcond, double *ferr, double *berr)
 Uses the factorization \(A = L D 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 tridiagonal matrix and X and B are n-by-nrhs matrices. More...
 
int64_t lapack::ptsvx (lapack::Factored fact, int64_t n, int64_t nrhs, double const *D, double const *E, double *DF, double *EF, double const *B, int64_t ldb, double *X, int64_t ldx, double *rcond, double *ferr, double *berr)
 
int64_t lapack::ptsvx (lapack::Factored fact, int64_t n, int64_t nrhs, float const *D, float const *E, float *DF, float *EF, float const *B, int64_t ldb, float *X, int64_t ldx, float *rcond, float *ferr, float *berr)
 
int64_t lapack::ptsvx (lapack::Factored fact, int64_t n, int64_t nrhs, float const *D, std::complex< float > const *E, float *DF, std::complex< float > *EF, std::complex< float > const *B, int64_t ldb, std::complex< float > *X, int64_t ldx, float *rcond, float *ferr, float *berr)
 

Detailed Description

Function Documentation

◆ ptsv()

int64_t lapack::ptsv ( int64_t  n,
int64_t  nrhs,
double *  D,
std::complex< double > *  E,
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 tridiagonal matrix, and X and B are n-by-nrhs matrices.

A is factored as \(A = L D L^H,\) and the factored form of A is then used to solve the system of equations.

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]nrhsThe number of right hand sides, i.e., the number of columns of the matrix B. nrhs >= 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 factorization \(A = L D L^H.\)
[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.
[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 is not positive definite, and the solution has not been computed. The factorization has not been completed unless i = n.

◆ ptsvx()

int64_t lapack::ptsvx ( lapack::Factored  fact,
int64_t  n,
int64_t  nrhs,
double const *  D,
std::complex< double > const *  E,
double *  DF,
std::complex< double > *  EF,
std::complex< double > const *  B,
int64_t  ldb,
std::complex< double > *  X,
int64_t  ldx,
double *  rcond,
double *  ferr,
double *  berr 
)

Uses the factorization \(A = L D 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 tridiagonal matrix 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]factSpecifies whether or not the factored form of the matrix A is supplied on entry.
  • lapack::Factored::Factored: On entry, DF and EF contain the factored form of A. D, E, DF, and EF will not be modified.
  • lapack::Factored::NotFactored: The matrix A will be copied to DF and EF and factored.
[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]DThe vector D of length n. The n diagonal elements of the tridiagonal matrix A.
[in]EThe vector E of length n-1. The (n-1) subdiagonal elements of the tridiagonal matrix A.
[in,out]DFThe vector DF of length n.
  • If fact = Factored, then DF is an input argument and on entry contains the n diagonal elements of the diagonal matrix D from the \(L D L^H\) factorization of A.
  • If fact = NotFactored, then DF is an output argument and on exit contains the n diagonal elements of the diagonal matrix D from the \(L D L^H\) factorization of A.
[in,out]EFThe vector EF of length n-1.
  • If fact = Factored, then EF is an input argument and on entry contains the (n-1) subdiagonal elements of the unit bidiagonal factor L from the \(L D L^H\) factorization of A.
  • If fact = NotFactored, then EF is an output argument and on exit contains the (n-1) subdiagonal elements of the unit bidiagonal factor L from the \(L D L^H\) factorization of A.
[in]BThe n-by-nrhs matrix B, stored in an ldb-by-nrhs array. The n-by-nrhs right hand side matrix 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.
[in]ldxThe leading dimension of the array X. ldx >= max(1,n).
[out]rcondThe reciprocal condition number of the matrix A. 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 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
> 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.