LAPACK
3.4.2
LAPACK: Linear Algebra PACKage
|
Functions/Subroutines | |
subroutine | cptcon (N, D, E, ANORM, RCOND, RWORK, INFO) |
CPTCON | |
subroutine | cpteqr (COMPZ, N, D, E, Z, LDZ, WORK, INFO) |
CPTEQR | |
subroutine | cptrfs (UPLO, N, NRHS, D, E, DF, EF, B, LDB, X, LDX, FERR, BERR, WORK, RWORK, INFO) |
CPTRFS | |
subroutine | cpttrf (N, D, E, INFO) |
CPTTRF | |
subroutine | cpttrs (UPLO, N, NRHS, D, E, B, LDB, INFO) |
CPTTRS | |
subroutine | cptts2 (IUPLO, N, NRHS, D, E, B, LDB) |
CPTTS2 solves a tridiagonal system of the form AX=B using the L D LH factorization computed by spttrf. |
This is the group of complex computational functions for PT matrices
subroutine cptcon | ( | integer | N, |
real, dimension( * ) | D, | ||
complex, dimension( * ) | E, | ||
real | ANORM, | ||
real | RCOND, | ||
real, dimension( * ) | RWORK, | ||
integer | INFO | ||
) |
CPTCON
Download CPTCON + dependencies [TGZ] [ZIP] [TXT]CPTCON computes the reciprocal of the condition number (in the 1-norm) of a complex Hermitian positive definite tridiagonal matrix using the factorization A = L*D*L**H or A = U**H*D*U computed by CPTTRF. 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 REAL array, dimension (N) The n diagonal elements of the diagonal matrix D from the factorization of A, as computed by CPTTRF. |
[in] | E | E is COMPLEX 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 CPTTRF. |
[in] | ANORM | ANORM is REAL The 1-norm of the original matrix A. |
[out] | RCOND | RCOND is REAL 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] | RWORK | RWORK is REAL 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 120 of file cptcon.f.
subroutine cpteqr | ( | character | COMPZ, |
integer | N, | ||
real, dimension( * ) | D, | ||
real, dimension( * ) | E, | ||
complex, dimension( ldz, * ) | Z, | ||
integer | LDZ, | ||
real, dimension( * ) | WORK, | ||
integer | INFO | ||
) |
CPTEQR
Download CPTEQR + dependencies [TGZ] [ZIP] [TXT]CPTEQR computes all eigenvalues and, optionally, eigenvectors of a symmetric positive definite tridiagonal matrix by first factoring the matrix using SPTTRF and then calling CBDSQR 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 positive definite Hermitian matrix can also be found if CHETRD, CHPTRD, or CHBTRD 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 Hermitian matrix also. Array Z contains the unitary 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 REAL 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 REAL 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 COMPLEX array, dimension (LDZ, N) On entry, if COMPZ = 'V', the unitary matrix used in the reduction to tridiagonal form. On exit, if COMPZ = 'V', the orthonormal eigenvectors of the original Hermitian 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 REAL 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 cpteqr.f.
subroutine cptrfs | ( | character | UPLO, |
integer | N, | ||
integer | NRHS, | ||
real, dimension( * ) | D, | ||
complex, dimension( * ) | E, | ||
real, dimension( * ) | DF, | ||
complex, dimension( * ) | EF, | ||
complex, dimension( ldb, * ) | B, | ||
integer | LDB, | ||
complex, dimension( ldx, * ) | X, | ||
integer | LDX, | ||
real, dimension( * ) | FERR, | ||
real, dimension( * ) | BERR, | ||
complex, dimension( * ) | WORK, | ||
real, dimension( * ) | RWORK, | ||
integer | INFO | ||
) |
CPTRFS
Download CPTRFS + dependencies [TGZ] [ZIP] [TXT]CPTRFS 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.
[in] | UPLO | UPLO is CHARACTER*1 Specifies whether the superdiagonal or the subdiagonal of the tridiagonal matrix A is stored and the form of the factorization: = 'U': E is the superdiagonal of A, and A = U**H*D*U; = 'L': E is the subdiagonal of A, and A = L*D*L**H. (The two forms are equivalent if A is real.) |
[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 REAL array, dimension (N) The n real diagonal elements of the tridiagonal matrix A. |
[in] | E | E is COMPLEX array, dimension (N-1) The (n-1) off-diagonal elements of the tridiagonal matrix A (see UPLO). |
[in] | DF | DF is REAL array, dimension (N) The n diagonal elements of the diagonal matrix D from the factorization computed by CPTTRF. |
[in] | EF | EF is COMPLEX array, dimension (N-1) The (n-1) off-diagonal elements of the unit bidiagonal factor U or L from the factorization computed by CPTTRF (see UPLO). |
[in] | B | B is COMPLEX 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 COMPLEX array, dimension (LDX,NRHS) On entry, the solution matrix X, as computed by CPTTRS. 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 REAL 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 REAL 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 COMPLEX array, dimension (N) |
[out] | RWORK | RWORK is REAL array, dimension (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 183 of file cptrfs.f.
subroutine cpttrf | ( | integer | N, |
real, dimension( * ) | D, | ||
complex, dimension( * ) | E, | ||
integer | INFO | ||
) |
CPTTRF
Download CPTTRF + dependencies [TGZ] [ZIP] [TXT]CPTTRF computes the L*D*L**H factorization of a complex Hermitian positive definite tridiagonal matrix A. The factorization may also be regarded as having the form A = U**H *D*U.
[in] | N | N is INTEGER The order of the matrix A. N >= 0. |
[in,out] | D | D is REAL 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**H factorization of A. |
[in,out] | E | E is COMPLEX 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**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. |
[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 93 of file cpttrf.f.
subroutine cpttrs | ( | character | UPLO, |
integer | N, | ||
integer | NRHS, | ||
real, dimension( * ) | D, | ||
complex, dimension( * ) | E, | ||
complex, dimension( ldb, * ) | B, | ||
integer | LDB, | ||
integer | INFO | ||
) |
CPTTRS
Download CPTTRS + dependencies [TGZ] [ZIP] [TXT]CPTTRS 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 CPTTRF. 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.
[in] | UPLO | UPLO is CHARACTER*1 Specifies 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. = 'U': A = U**H*D*U, E is the superdiagonal of U = 'L': A = L*D*L**H, E is the subdiagonal of L |
[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 REAL array, dimension (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] | E | E is COMPLEX array, dimension (N-1) If UPLO = 'U', the (n-1) superdiagonal elements of the unit bidiagonal factor U from the factorization A = U**H*D*U. If UPLO = 'L', the (n-1) subdiagonal elements of the unit bidiagonal factor L from the factorization A = L*D*L**H. |
[in,out] | B | B is REAL 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 122 of file cpttrs.f.
subroutine cptts2 | ( | integer | IUPLO, |
integer | N, | ||
integer | NRHS, | ||
real, dimension( * ) | D, | ||
complex, dimension( * ) | E, | ||
complex, dimension( ldb, * ) | B, | ||
integer | LDB | ||
) |
CPTTS2 solves a tridiagonal system of the form AX=B using the L D LH factorization computed by spttrf.
Download CPTTS2 + dependencies [TGZ] [ZIP] [TXT]CPTTS2 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 CPTTRF. 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.
[in] | IUPLO | IUPLO is INTEGER Specifies 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. = 1: A = U**H *D*U, E is the superdiagonal of U = 0: A = L*D*L**H, E is the subdiagonal of L |
[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 REAL array, dimension (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] | E | E is COMPLEX array, dimension (N-1) If IUPLO = 1, the (n-1) superdiagonal elements of the unit bidiagonal factor U from the factorization A = U**H*D*U. If IUPLO = 0, the (n-1) subdiagonal elements of the unit bidiagonal factor L from the factorization A = L*D*L**H. |
[in,out] | B | B is REAL 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 114 of file cptts2.f.