F07HAF (DPBSV) (PDF version)
F07 Chapter Contents
F07 Chapter Introduction
NAG Library Manual

NAG Library Routine Document

F07HAF (DPBSV)

Note:  before using this routine, please read the Users' Note for your implementation to check the interpretation of bold italicised terms and other implementation-dependent details.

+ Contents

    1  Purpose
    7  Accuracy

1  Purpose

F07HAF (DPBSV) computes the solution to a real system of linear equations
AX=B ,
where A is an n by n symmetric positive definite band matrix of bandwidth 2 kd + 1  and X and B are n by r matrices.

2  Specification

SUBROUTINE F07HAF ( UPLO, N, KD, NRHS, AB, LDAB, B, LDB, INFO)
INTEGER  N, KD, NRHS, LDAB, LDB, INFO
REAL (KIND=nag_wp)  AB(LDAB,*), B(LDB,*)
CHARACTER(1)  UPLO
The routine may be called by its LAPACK name dpbsv.

3  Description

F07HAF (DPBSV) uses the Cholesky decomposition to factor A as A=UTU if UPLO='U' or A=LLT if UPLO='L', where U is an upper triangular band matrix, and L is a lower triangular band matrix, with the same number of superdiagonals or subdiagonals as A. The factored form of A is then used to solve the system of equations AX=B.

4  References

Anderson E, Bai Z, Bischof C, Blackford S, Demmel J, Dongarra J J, Du Croz J J, Greenbaum A, Hammarling S, McKenney A and Sorensen D (1999) LAPACK Users' Guide (3rd Edition) SIAM, Philadelphia http://www.netlib.org/lapack/lug
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore

5  Parameters

1:     UPLO – CHARACTER(1)Input
On entry: if UPLO='U', the upper triangle of A is stored.
If UPLO='L', the lower triangle of A is stored.
Constraint: UPLO='U' or 'L'.
2:     N – INTEGERInput
On entry: n, the number of linear equations, i.e., the order of the matrix A.
Constraint: N0.
3:     KD – INTEGERInput
On entry: kd, the number of superdiagonals of the matrix A if UPLO='U', or the number of subdiagonals if UPLO='L'.
Constraint: KD0.
4:     NRHS – INTEGERInput
On entry: r, the number of right-hand sides, i.e., the number of columns of the matrix B.
Constraint: NRHS0.
5:     AB(LDAB,*) – REAL (KIND=nag_wp) arrayInput/Output
Note: the second dimension of the array AB must be at least max1,N.
On entry: the upper or lower triangle of the symmetric band matrix A.
The matrix is stored in rows 1 to kd+1, more precisely,
  • if UPLO='U', the elements of the upper triangle of A within the band must be stored with element Aij in ABkd+1+i-jj​ for ​max1,j-kdij;
  • if UPLO='L', the elements of the lower triangle of A within the band must be stored with element Aij in AB1+i-jj​ for ​jiminn,j+kd.
On exit: if INFO=0, the triangular factor U or L from the Cholesky factorization A=UTU or A=LLT of the band matrix A, in the same storage format as A.
6:     LDAB – INTEGERInput
On entry: the first dimension of the array AB as declared in the (sub)program from which F07HAF (DPBSV) is called.
Constraint: LDABKD+1.
7:     B(LDB,*) – REAL (KIND=nag_wp) arrayInput/Output
Note: the second dimension of the array B must be at least max1,NRHS.
On entry: the n by r right-hand side matrix B.
On exit: if INFO=0, the n by r solution matrix X.
8:     LDB – INTEGERInput
On entry: the first dimension of the array B as declared in the (sub)program from which F07HAF (DPBSV) is called.
Constraint: LDBmax1,N.
9:     INFO – INTEGEROutput
On exit: INFO=0 unless the routine detects an error (see Section 6).

6  Error Indicators and Warnings

Errors or warnings detected by the routine:
INFO<0
If INFO=-i, the ith argument had an illegal value. An explanatory message is output, and execution of the program is terminated.
INFO>0
If INFO=i, 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.

7  Accuracy

The computed solution for a single right-hand side, x^ , satisfies an equation of the form
A+E x^=b ,
where
E1 = Oε A1
and ε  is the machine precision. An approximate error bound for the computed solution is given by
x^-x1 x1 κA E1 A1 ,
where κA = A-11 A1 , the condition number of A  with respect to the solution of the linear equations. See Section 4.4 of Anderson et al. (1999) for further details.
F07HBF (DPBSVX) is a comprehensive LAPACK driver that returns forward and backward error bounds and an estimate of the condition number. Alternatively, F04BFF solves Ax=b  and returns a forward error bound and condition estimate. F04BFF calls F07HAF (DPBSV) to solve the equations.

8  Further Comments

When nk , the total number of floating point operations is approximately nk+12+4nkr , where k  is the number of superdiagonals and r  is the number of right-hand sides.
The complex analogue of this routine is F07HNF (ZPBSV).

9  Example

This example solves the equations
Ax=b ,
where A  is the symmetric positive definite band matrix
A = 5.49 2.68 0.00 0.00 2.68 5.63 -2.39 0.00 0.00 -2.39 2.60 -2.22 0.00 0.00 -2.22 5.17   and   b = 22.09 9.31 -5.24 11.83 .
Details of the Cholesky factorization of A  are also output.

9.1  Program Text

Program Text (f07hafe.f90)

9.2  Program Data

Program Data (f07hafe.d)

9.3  Program Results

Program Results (f07hafe.r)


F07HAF (DPBSV) (PDF version)
F07 Chapter Contents
F07 Chapter Introduction
NAG Library Manual

© The Numerical Algorithms Group Ltd, Oxford, UK. 2012