|
|
◆ dposv()
| def dposv |
( |
uplo |
, |
|
|
n |
, |
|
|
a |
, |
|
|
b |
, |
|
|
nrhs |
= 1 |
|
) |
| |
Solution to system of linear equations AX = B for a symmetric positive definite matrix
- Purpose
- dposv computes the solution to a real system of linear equations where A is an n x n symmetric positive definite matrix and X and B are n x nrhs matrices.
The Cholesky decomposition is used to factor A as A = U^T*U, if uplo = 'U', or
A = L*L^T, if uplo = 'L',
where U is an upper triangular matrix and L is a lower triangular matrix. The factored form of A is then used to solve the system of equations A * X = B.
- Returns
- info (int)
= 0: Successful exit
= -1: The argument uplo had an illegal value (uplo != 'U' nor 'L')
= -2: The argument n had an illegal value (n < 0)
= -3: The argument a is invalid.
= -4: The argument b is invalid.
= -5: The argument nrhs had an illegal value (nrhs < 0)
= i > 0: 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.
- Parameters
-
| [in] | uplo | = 'U': Upper triangle of A is stored.
= 'L': Lower triangle of A is stored. |
| [in] | n | Number of linear equations, i.e., order of the matrix A. (n >= 0) (If n = 0, returns without computation) |
| [in,out] | a | Numpy ndarray (2-dimensional, float, n x n)
[in] n x n symmetric positove definite matrix A. The upper or lower triangular part is to be referenced in accordance with uplo.
[out] If info = 0, the factor U or L from the Cholesky factorization A = U^T*U or A = L*L^T. |
| [in,out] | b | Numpy ndarray (1 or 2-dimensional, float, n or n x nrhs)
[in] n x nrhs right hand side matrix B.
[out] If info = 0, the n x nrhs solution matrix X. |
| [in] | nrhs | (Optional)
Number of right hand sides, i.e., number of columns of the matrix B. (nrhs >= 0) (If nrhs = 0, returns without computation) (default = 1) |
- Reference
- LAPACK
- Example Program
- Solve the system of linear equations Ax = B and estimate the reciprocal of the condition number (RCond) of A, where A is symmetric positive definite and
( 2.2 -0.11 -0.32 ) ( -1.566 )
A = ( -0.11 2.93 0.81 ), B = ( -2.8425 )
( -0.32 0.81 2.37 ) ( -1.1765 )
def TestDposv():
n = 3
a = np.array([
[2.2, 0.0, 0.0],
[-0.11, 2.93, 0.0],
[-0.32, 0.81, 2.37]])
b = np.array([-1.566, -2.8425, -1.1765])
anorm, info = dlansy( '1', 'U', n, a)
info = dposv( 'U', n, a, b)
print(b, info)
rcond, info = dpocon( 'U', n, a, anorm)
print(rcond, info)
def dposv(uplo, n, a, b, nrhs=1) Solution to system of linear equations AX = B for a symmetric positive definite matrix
def dpocon(uplo, n, a, anorm) Condition number of a symmetric positive definite matrix
def dlansy(norm, uplo, n, a) One norm, Frobenius norm, infinity norm, or largest absolute value of any element of a real symmetric...
- Example Results
>>> TestDposv()
[-0.8 -0.92 -0.29] 0
0.4467910780689557 0
|