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

◆ Zpbsv()

Sub Zpbsv ( Uplo As  String,
N As  Long,
Kd As  Long,
Ab() 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 band matrix

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 band matrix, 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 band matrix, and L is a lower triangular band matrix, with the same number of super-diagonals or sub-diagonals as A. 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]KdNumber of super-diagonals of the matrix A if Uplo = "U", or number of sub-diagonals if Uplo = "L". (Kd >= 0)
[in,out]Ab()Array Ab(LAb1 - 1, LAb2 - 1) (LAb1 >= Kd + 1, LAb2 >= N)
[in] N x N Hermitian positive definite band matrix A in Kd+1 x N symmetric band matrix form. Upper or lower part is to be stored in accordance with uplo. See below for further details.
[out] If Info = 0, the triangular factor U or L from the Cholesky factorization A = U^H*U or A = L*L^H of the band matrix A, in the same storage format as A.
[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 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 Kd had an illegal value. (Kd < 0)
= -4: The argument Ab() is invalid.
= -5: The argument B() is invalid.
= -7: 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)
Further Details
The symmetric band matrix form is illustrated by the following example, when n = 6, kd = 2, and uplo = "U":
On entry:

   *    *   a13  a24  a35  a46
   *   a12  a23  a34  a45  a56
  a11  a22  a33  a44  a55  a66

On exit:

   *    *   u13  u24  u35  u46
   *   u12  u23  u34  u45  u56
  u11  u22  u33  u44  u55  u66
Similarly, if uplo = "L" the format of A is as follows:
On entry:

  a11  a22  a33  a44  a55  a66
  a21  a32  a43  a54  a65   *
  a31  a42  a53  a64   *    *

On exit:

  l11  l22  l33  l44  l55  l66
  l21  l32  l43  l54  l65   *
  l31  l42  l53  l64   *    *
Array elements marked * are not used by the routine.
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.88 0.29-0.44i 0 )
A = ( 0.29+0.44i 0.62 -0.01-0.02i )
( 0 -0.01+0.02i 0.46 )
( 1.6236-0.7300i )
B = ( 0.1581+0.1537i )
( 0.1132-0.2290i )
Sub Ex_Zpbsv()
Const N As Long = 3, Kd = 1
Dim Ab(Kd, N - 1) As Complex, B(N - 1) As Complex
Dim ANorm As Double, RCond As Double, Info As Long
Ab(0, 0) = cmplx(2.88): Ab(0, 1) = cmplx(0.62): Ab(0, 2) = cmplx(0.46)
Ab(1, 0) = cmplx(0.29, 0.44): Ab(1, 1) = cmplx(-0.01, 0.02)
B(0) = cmplx(1.6236, -0.73): B(1) = cmplx(0.1581, 0.1537): B(2) = cmplx(0.1132, -0.229)
ANorm = Zlanhb("1", "L", N, Kd, Ab())
Call Zpbsv("L", N, Kd, Ab(), B(), Info)
If Info = 0 Then Call Zpbcon("L", N, Kd, Ab(), 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 Cimag(A As Complex) As Double
Imaginary part of complex number
Function Creal(A As Complex) As Double
Real part of complex number
Function Zlanhb(Norm As String, Uplo As String, N As Long, K As Long, Ab() As Complex, Optional Info As Long) As Double
One norm, Frobenius norm, infinity norm, or largest absolute value of any element of a Hermitian band...
Sub Zpbcon(Uplo As String, N As Long, Kd As Long, Ab() As Complex, ANorm As Double, RCond As Double, Info As Long)
Condition number of a Hermitian positive definite band matrix
Sub Zpbsv(Uplo As String, N As Long, Kd As Long, Ab() 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 band ...
Example Results
X = 0.59 -0.28 -0.2 -0.04 0.24 -0.49
RCond = 0.124521368143895
Info = 0