LAPACK
3.4.2
LAPACK: Linear Algebra PACKage
|
Functions/Subroutines | |
subroutine | dptcon (N, D, E, ANORM, RCOND, WORK, INFO) |
DPTCON | |
subroutine | dpteqr (COMPZ, N, D, E, Z, LDZ, WORK, INFO) |
DPTEQR | |
subroutine | dptrfs (N, NRHS, D, E, DF, EF, B, LDB, X, LDX, FERR, BERR, WORK, INFO) |
DPTRFS | |
subroutine | dpttrf (N, D, E, INFO) |
DPTTRF | |
subroutine | dpttrs (N, NRHS, D, E, B, LDB, INFO) |
DPTTRS | |
subroutine | dptts2 (N, NRHS, D, E, B, LDB) |
DPTTS2 solves a tridiagonal system of the form AX=B using the L D LH factorization computed by spttrf. |
This is the group of double computational functions for PT matrices
subroutine dptcon | ( | integer | N, |
double precision, dimension( * ) | D, | ||
double precision, dimension( * ) | E, | ||
double precision | ANORM, | ||
double precision | RCOND, | ||
double precision, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
DPTCON
Download DPTCON + dependencies [TGZ] [ZIP] [TXT]DPTCON computes the reciprocal of the condition number (in the 1-norm) of a real symmetric positive definite tridiagonal matrix using the factorization A = L*D*L**T or A = U**T*D*U computed by DPTTRF. Norm(inv(A)) is computed by a direct method, and the reciprocal of the condition number is computed as RCOND = 1 / (ANORM * norm(inv(A))).
[in] | N | N is INTEGER The order of the matrix A. N >= 0. |
[in] | D | D is DOUBLE PRECISION array, dimension (N) The n diagonal elements of the diagonal matrix D from the factorization of A, as computed by DPTTRF. |
[in] | E | E is DOUBLE PRECISION array, dimension (N-1) The (n-1) off-diagonal elements of the unit bidiagonal factor U or L from the factorization of A, as computed by DPTTRF. |
[in] | ANORM | ANORM is DOUBLE PRECISION The 1-norm of the original matrix A. |
[out] | RCOND | RCOND is DOUBLE PRECISION The reciprocal of the condition number of the matrix A, computed as RCOND = 1/(ANORM * AINVNM), where AINVNM is the 1-norm of inv(A) computed in this routine. |
[out] | WORK | WORK is DOUBLE PRECISION array, dimension (N) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
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.
Definition at line 119 of file dptcon.f.
subroutine dpteqr | ( | character | COMPZ, |
integer | N, | ||
double precision, dimension( * ) | D, | ||
double precision, dimension( * ) | E, | ||
double precision, dimension( ldz, * ) | Z, | ||
integer | LDZ, | ||
double precision, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
DPTEQR
Download DPTEQR + dependencies [TGZ] [ZIP] [TXT]DPTEQR computes all eigenvalues and, optionally, eigenvectors of a symmetric positive definite tridiagonal matrix by first factoring the matrix using DPTTRF, and then calling DBDSQR to compute the singular values of the bidiagonal factor. This routine computes the eigenvalues of the positive definite tridiagonal matrix to high relative accuracy. This means that if the eigenvalues range over many orders of magnitude in size, then the small eigenvalues and corresponding eigenvectors will be computed more accurately than, for example, with the standard QR method. The eigenvectors of a full or band symmetric positive definite matrix can also be found if DSYTRD, DSPTRD, or DSBTRD has been used to reduce this matrix to tridiagonal form. (The reduction to tridiagonal form, however, may preclude the possibility of obtaining high relative accuracy in the small eigenvalues of the original matrix, if these eigenvalues range over many orders of magnitude.)
[in] | COMPZ | COMPZ is CHARACTER*1 = 'N': Compute eigenvalues only. = 'V': Compute eigenvectors of original symmetric matrix also. Array Z contains the orthogonal matrix used to reduce the original matrix to tridiagonal form. = 'I': Compute eigenvectors of tridiagonal matrix also. |
[in] | N | N is INTEGER The order of the matrix. N >= 0. |
[in,out] | D | D is DOUBLE PRECISION array, dimension (N) On entry, the n diagonal elements of the tridiagonal matrix. On normal exit, D contains the eigenvalues, in descending order. |
[in,out] | E | E is DOUBLE PRECISION array, dimension (N-1) On entry, the (n-1) subdiagonal elements of the tridiagonal matrix. On exit, E has been destroyed. |
[in,out] | Z | Z is DOUBLE PRECISION array, dimension (LDZ, N) On entry, if COMPZ = 'V', the orthogonal matrix used in the reduction to tridiagonal form. On exit, if COMPZ = 'V', the orthonormal eigenvectors of the original symmetric matrix; if COMPZ = 'I', the orthonormal eigenvectors of the tridiagonal matrix. If INFO > 0 on exit, Z contains the eigenvectors associated with only the stored eigenvalues. If COMPZ = 'N', then Z is not referenced. |
[in] | LDZ | LDZ is INTEGER The leading dimension of the array Z. LDZ >= 1, and if COMPZ = 'V' or 'I', LDZ >= max(1,N). |
[out] | WORK | WORK is DOUBLE PRECISION array, dimension (4*N) |
[out] | INFO | INFO is INTEGER = 0: successful exit. < 0: if INFO = -i, the i-th argument had an illegal value. > 0: if INFO = i, and i is: <= N the Cholesky factorization of the matrix could not be performed because the i-th principal minor was not positive definite. > N the SVD algorithm failed to converge; if INFO = N+i, i off-diagonal elements of the bidiagonal factor did not converge to zero. |
Definition at line 146 of file dpteqr.f.
subroutine dptrfs | ( | integer | N, |
integer | NRHS, | ||
double precision, dimension( * ) | D, | ||
double precision, dimension( * ) | E, | ||
double precision, dimension( * ) | DF, | ||
double precision, dimension( * ) | EF, | ||
double precision, dimension( ldb, * ) | B, | ||
integer | LDB, | ||
double precision, dimension( ldx, * ) | X, | ||
integer | LDX, | ||
double precision, dimension( * ) | FERR, | ||
double precision, dimension( * ) | BERR, | ||
double precision, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
DPTRFS
Download DPTRFS + dependencies [TGZ] [ZIP] [TXT]DPTRFS improves the computed solution to a system of linear equations when the coefficient matrix is symmetric positive definite and tridiagonal, and provides error bounds and backward error estimates for the solution.
[in] | N | N is INTEGER The order of the matrix A. N >= 0. |
[in] | NRHS | NRHS is INTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0. |
[in] | D | D is DOUBLE PRECISION array, dimension (N) The n diagonal elements of the tridiagonal matrix A. |
[in] | E | E is DOUBLE PRECISION array, dimension (N-1) The (n-1) subdiagonal elements of the tridiagonal matrix A. |
[in] | DF | DF is DOUBLE PRECISION array, dimension (N) The n diagonal elements of the diagonal matrix D from the factorization computed by DPTTRF. |
[in] | EF | EF is DOUBLE PRECISION array, dimension (N-1) The (n-1) subdiagonal elements of the unit bidiagonal factor L from the factorization computed by DPTTRF. |
[in] | B | B is DOUBLE PRECISION array, dimension (LDB,NRHS) The right hand side matrix B. |
[in] | LDB | LDB is INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[in,out] | X | X is DOUBLE PRECISION array, dimension (LDX,NRHS) On entry, the solution matrix X, as computed by DPTTRS. On exit, the improved solution matrix X. |
[in] | LDX | LDX is INTEGER The leading dimension of the array X. LDX >= max(1,N). |
[out] | FERR | FERR is DOUBLE PRECISION array, dimension (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] | BERR | BERR is DOUBLE PRECISION array, dimension (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). |
[out] | WORK | WORK is DOUBLE PRECISION array, dimension (2*N) |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -i, the i-th argument had an illegal value |
ITMAX is the maximum number of steps of iterative refinement.
Definition at line 163 of file dptrfs.f.
subroutine dpttrf | ( | integer | N, |
double precision, dimension( * ) | D, | ||
double precision, dimension( * ) | E, | ||
integer | INFO | ||
) |
DPTTRF
Download DPTTRF + dependencies [TGZ] [ZIP] [TXT]DPTTRF computes the L*D*L**T factorization of a real symmetric positive definite tridiagonal matrix A. The factorization may also be regarded as having the form A = U**T*D*U.
[in] | N | N is INTEGER The order of the matrix A. N >= 0. |
[in,out] | D | D is DOUBLE PRECISION array, dimension (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**T factorization of A. |
[in,out] | E | E is DOUBLE PRECISION array, dimension (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**T factorization of A. E can also be regarded as the superdiagonal of the unit bidiagonal factor U from the U**T*D*U factorization of A. |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -k, the k-th argument had an illegal value > 0: if INFO = k, the leading minor of order k is not positive definite; if k < N, the factorization could not be completed, while if k = N, the factorization was completed, but D(N) <= 0. |
Definition at line 92 of file dpttrf.f.
subroutine dpttrs | ( | integer | N, |
integer | NRHS, | ||
double precision, dimension( * ) | D, | ||
double precision, dimension( * ) | E, | ||
double precision, dimension( ldb, * ) | B, | ||
integer | LDB, | ||
integer | INFO | ||
) |
DPTTRS
Download DPTTRS + dependencies [TGZ] [ZIP] [TXT]DPTTRS solves a tridiagonal system of the form A * X = B using the L*D*L**T factorization of A computed by DPTTRF. D is a diagonal matrix specified in the vector D, L is a unit bidiagonal matrix whose subdiagonal is specified in the vector E, and X and B are N by NRHS matrices.
[in] | N | N is INTEGER The order of the tridiagonal matrix A. N >= 0. |
[in] | NRHS | NRHS is INTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0. |
[in] | D | D is DOUBLE PRECISION array, dimension (N) The n diagonal elements of the diagonal matrix D from the L*D*L**T factorization of A. |
[in] | E | E is DOUBLE PRECISION array, dimension (N-1) The (n-1) subdiagonal elements of the unit bidiagonal factor L from the L*D*L**T factorization of A. E can also be regarded as the superdiagonal of the unit bidiagonal factor U from the factorization A = U**T*D*U. |
[in,out] | B | B is DOUBLE PRECISION array, dimension (LDB,NRHS) On entry, the right hand side vectors B for the system of linear equations. On exit, the solution vectors, X. |
[in] | LDB | LDB is INTEGER The leading dimension of the array B. LDB >= max(1,N). |
[out] | INFO | INFO is INTEGER = 0: successful exit < 0: if INFO = -k, the k-th argument had an illegal value |
Definition at line 110 of file dpttrs.f.
subroutine dptts2 | ( | integer | N, |
integer | NRHS, | ||
double precision, dimension( * ) | D, | ||
double precision, dimension( * ) | E, | ||
double precision, dimension( ldb, * ) | B, | ||
integer | LDB | ||
) |
DPTTS2 solves a tridiagonal system of the form AX=B using the L D LH factorization computed by spttrf.
Download DPTTS2 + dependencies [TGZ] [ZIP] [TXT]DPTTS2 solves a tridiagonal system of the form A * X = B using the L*D*L**T factorization of A computed by DPTTRF. D is a diagonal matrix specified in the vector D, L is a unit bidiagonal matrix whose subdiagonal is specified in the vector E, and X and B are N by NRHS matrices.
[in] | N | N is INTEGER The order of the tridiagonal matrix A. N >= 0. |
[in] | NRHS | NRHS is INTEGER The number of right hand sides, i.e., the number of columns of the matrix B. NRHS >= 0. |
[in] | D | D is DOUBLE PRECISION array, dimension (N) The n diagonal elements of the diagonal matrix D from the L*D*L**T factorization of A. |
[in] | E | E is DOUBLE PRECISION array, dimension (N-1) The (n-1) subdiagonal elements of the unit bidiagonal factor L from the L*D*L**T factorization of A. E can also be regarded as the superdiagonal of the unit bidiagonal factor U from the factorization A = U**T*D*U. |
[in,out] | B | B is DOUBLE PRECISION array, dimension (LDB,NRHS) On entry, the right hand side vectors B for the system of linear equations. On exit, the solution vectors, X. |
[in] | LDB | LDB is INTEGER The leading dimension of the array B. LDB >= max(1,N). |
Definition at line 103 of file dptts2.f.