F07HDF (DPBTRF) (PDF version)
F07 Chapter Contents
F07 Chapter Introduction
NAG Library Manual

NAG Library Routine Document

F07HDF (DPBTRF)

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

F07HDF (DPBTRF) computes the Cholesky factorization of a real symmetric positive definite band matrix.

2  Specification

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

3  Description

F07HDF (DPBTRF) forms the Cholesky factorization of a real symmetric positive definite band matrix A either as A=UTU if UPLO='U' or A=LLT if UPLO='L', where U (or L) is an upper (or lower) triangular band matrix with the same number of superdiagonals (or subdiagonals) as A.

4  References

Demmel J W (1989) On floating-point errors in Cholesky LAPACK Working Note No. 14 University of Tennessee, Knoxville
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: specifies whether the upper or lower triangular part of A is stored and how A is to be factorized.
UPLO='U'
The upper triangular part of A is stored and A is factorized as UTU, where U is upper triangular.
UPLO='L'
The lower triangular part of A is stored and A is factorized as LLT, where L is lower triangular.
Constraint: UPLO='U' or 'L'.
2:     N – INTEGERInput
On entry: n, the order of the matrix A.
Constraint: N0.
3:     KD – INTEGERInput
On entry: kd, the number of superdiagonals or subdiagonals of the matrix A.
Constraint: KD0.
4:     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 n by n 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: the upper or lower triangle of A is overwritten by the Cholesky factor U or L as specified by UPLO, using the same storage format as described above.
5:     LDAB – INTEGERInput
On entry: the first dimension of the array AB as declared in the (sub)program from which F07HDF (DPBTRF) is called.
Constraint: LDABKD+1.
6:     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 parameter 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 is not positive definite and the factorization could not be completed. Hence A itself is not positive definite. This may indicate an error in forming the matrix A. There is no routine specifically designed to factorize a band matrix which is not positive definite; the matrix must be treated either as a nonsymmetric band matrix, by calling F07BDF (DGBTRF) or as a full matrix, by calling F07MDF (DSYTRF).

7  Accuracy

If UPLO='U', the computed factor U is the exact factor of a perturbed matrix A+E, where
Eck+1εUTU ,
ck+1 is a modest linear function of k+1, and ε is the machine precision.
If UPLO='L', a similar statement holds for the computed factor L. It follows that eijck+1εaiiajj.

8  Further Comments

The total number of floating point operations is approximately n k+1 2, assuming nk.
A call to F07HDF (DPBTRF) may be followed by calls to the routines:
The complex analogue of this routine is F07HRF (ZPBTRF).

9  Example

This example computes the Cholesky factorization of the matrix A, where
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 .

9.1  Program Text

Program Text (f07hdfe.f90)

9.2  Program Data

Program Data (f07hdfe.d)

9.3  Program Results

Program Results (f07hdfe.r)


F07HDF (DPBTRF) (PDF version)
F07 Chapter Contents
F07 Chapter Introduction
NAG Library Manual

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