XLPack 7.0
XLPack Numerical Library (Excel VBA) Reference Manual
Loading...
Searching...
No Matches

◆ Zppsv()

Sub Zppsv ( Uplo As  String,
N As  Long,
Ap() As  Complex,
B() As  Complex,
Info As  Long,
Optional Nrhs As  Long = 1 
)

(Simple driver) Solution to system of linear equations AX = B for a Hermitian positive definite matrix in packed form

Purpose
This routine computes the solution to a complex system of linear equations
A * X = B,
where A is an n x n Hermitian positive definite matrix stored in packed form and X and B are n x nrhs matrices.

The Cholesky decomposition is used to factor A as
A = U^H * U, if Uplo = "U", or
A = L * L^H, 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.
Parameters
[in]Uplo= "U": Upper triangle of A is stored.
= "L": Lower triangle of A is stored.
[in]NNumber of linear equations, i.e., order of the matrix A. (N >= 0) (If N = 0, returns without computation)
[in,out]Ap()Array Ap(LAp - 1) (LAp >= N(N + 1)/2)
[in] N x N Hermitian positive definite matrix A in packed form. The upper or lower part is to be stored in accordance with Uplo.
[out] If Info = 0, the factor U or L from the Cholesky factorization A = U^H*U or A = L*L^H.
[in,out]B()Array B(LB1 - 1, LB2 - 1) (LB1 >= max(1, N), LB2 >= Nrhs) (2D array) or B(LB - 1) (LB >= max(1, N), Nrhs = 1) (1D array)
[in] N x Nrhs matrix of right hand side matrix B.
[out] If Info = 0, the N x Nrhs solution matrix X.
[out]Info= 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 Ap() is invalid.
= -4: The argument B() is invalid.
= -6: 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.
[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
( 2.20 -0.11+0.93i 0.81-0.37i )
A = ( -0.11-0.93i 2.32 -0.80+0.92i )
( 0.81+0.37i -0.80-0.92i 2.29 )
( 1.5980+1.4644i )
B = ( 1.3498+1.4398i )
( 2.0561-0.5441i )
Sub Ex_Zppsv()
Const N As Long = 3
Dim Ap(N * (N + 1) / 2) As Complex, B(N - 1) As Complex
Dim ANorm As Double, RCond As Double, Info As Long
Ap(0) = Cmplx(2.2, 0)
Ap(1) = Cmplx(-0.11, -0.93): Ap(3) = Cmplx(2.32, 0)
Ap(2) = Cmplx(0.81, 0.37): Ap(4) = Cmplx(-0.8, -0.92): Ap(5) = Cmplx(2.29, 0)
B(0) = Cmplx(1.598, 1.4644): B(1) = Cmplx(1.3498, 1.4398): B(2) = Cmplx(2.0561, -0.5441)
ANorm = Zlanhp("1", "L", N, Ap())
Call Zppsv("L", N, Ap(), B(), Info)
If Info = 0 Then Call Zppcon("L", N, Ap(), ANorm, RCond, Info)
Debug.Print "X =",
Debug.Print Creal(B(0)), Cimag(B(0)), Creal(B(1)), Cimag(B(1)), Creal(B(2)), Cimag(B(2))
Debug.Print "RCond =", RCond
Debug.Print "Info =", Info
End Sub
Function Cmplx(R As Double, Optional I As Double=0) As Complex
Building complex number
Function Cimag(A As Complex) As Double
Imaginary part of complex number
Function Creal(A As Complex) As Double
Real part of complex number
Function Zlanhp(Norm As String, Uplo As String, N As Long, Ap() As Complex, Optional Info As Long) As Double
One norm, Frobenius norm, infinity norm, or largest absolute value of any element of a Hermitian matr...
Sub Zppcon(Uplo As String, N As Long, Ap() As Complex, ANorm As Double, RCond As Double, Info As Long)
Condition number of a Hermitian positive definite matrix in packed form
Sub Zppsv(Uplo As String, N As Long, Ap() As Complex, B() As Complex, Info As Long, Optional Nrhs As Long=1)
(Simple driver) Solution to system of linear equations AX = B for a Hermitian positive definite matri...
Example Results
X = 0.86 0.64 0.51 0.71 0.59 -0.15
RCond = 8.85790434328042E-02
Info = 0