hide long namesshow long names
hide short namesshow short names
Integer type:  int32  int64  nag_int  show int32  show int32  show int64  show int64  show nag_int  show nag_int

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

NAG Toolbox: nag_lapack_dppsv (f07ga)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_lapack_dppsv (f07ga) computes the solution to a real system of linear equations
AX=B ,  
where A is an n by n symmetric positive definite matrix stored in packed format and X and B are n by r matrices.

Syntax

[ap, b, info] = f07ga(uplo, ap, b, 'n', n, 'nrhs_p', nrhs_p)
[ap, b, info] = nag_lapack_dppsv(uplo, ap, b, 'n', n, 'nrhs_p', nrhs_p)

Description

nag_lapack_dppsv (f07ga) 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 matrix and L is a lower triangular matrix. The factored form of A is then used to solve the system of equations AX=B.

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

Parameters

Compulsory Input Parameters

1:     uplo – string (length ≥ 1)
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:     ap: – double array
The dimension of the array ap must be at least max1,n×n+1/2
The n by n symmetric matrix A, packed by columns.
More precisely,
  • if uplo='U', the upper triangle of A must be stored with element Aij in api+jj-1/2 for ij;
  • if uplo='L', the lower triangle of A must be stored with element Aij in api+2n-jj-1/2 for ij.
3:     bldb: – double array
The first dimension of the array b must be at least max1,n.
The second dimension of the array b must be at least max1,nrhs_p.
The n by r right-hand side matrix B.

Optional Input Parameters

1:     n int64int32nag_int scalar
Default: the first dimension of the array b.
n, the number of linear equations, i.e., the order of the matrix A.
Constraint: n0.
2:     nrhs_p int64int32nag_int scalar
Default: the second dimension of the array b.
r, the number of right-hand sides, i.e., the number of columns of the matrix B.
Constraint: nrhs_p0.

Output Parameters

1:     ap: – double array
The dimension of the array ap will be max1,n×n+1/2
If info=0, the factor U or L from the Cholesky factorization A=UTU or A=LLT, in the same storage format as A.
2:     bldb: – double array
The first dimension of the array b will be max1,n.
The second dimension of the array b will be max1,nrhs_p.
If info=0, the n by r solution matrix X.
3:     info int64int32nag_int scalar
info=0 unless the function detects an error (see Error Indicators and Warnings).

Error Indicators and Warnings

   info<0
If info=-i, argument i had an illegal value. An explanatory message is output, and execution of the program is terminated.
   info>0
The leading minor of order _ of A is not positive definite, so the factorization could not be completed, and the solution has not been computed.

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.
nag_lapack_dppsvx (f07gb) is a comprehensive LAPACK driver that returns forward and backward error bounds and an estimate of the condition number. Alternatively, nag_linsys_real_posdef_packed_solve (f04be) solves Ax=b  and returns a forward error bound and condition estimate. nag_linsys_real_posdef_packed_solve (f04be) calls nag_lapack_dppsv (f07ga) to solve the equations.

Further Comments

The total number of floating-point operations is approximately 13 n3 + 2n2r , where r  is the number of right-hand sides.
The complex analogue of this function is nag_lapack_zppsv (f07gn).

Example

This example solves the equations
Ax=b ,  
where A  is the symmetric positive definite matrix
A = 4.16 -3.12 0.56 -0.10 -3.12 5.03 -0.83 1.18 0.56 -0.83 0.76 0.34 -0.10 1.18 0.34 1.18   and   b = 8.70 -13.35 1.89 -4.14 .  
Details of the Cholesky factorization of A  are also output.
function f07ga_example


fprintf('f07ga example results\n\n');

% Symmetric matrix A, upper triangular part packed in ap
uplo = 'U';
ap = [4.16  ...
     -3.12  5.03  ...
      0.56 -0.83  0.76  ...
     -0.10  1.18  0.34  1.18];

% Rhs
b = [  8.7;
     -13.35;
      1.89;
     -4.14];
n = int64(size(b,1));

% Solve
[apf, x, info] = f07ga( ...
                        uplo, ap, b);

disp('Solution');
disp(x');

[ifail] = x04cc( ...
                 uplo, 'N', n, apf, 'Cholesky factor');


f07ga example results

Solution
    1.0000   -1.0000    2.0000   -3.0000

 Cholesky factor
             1          2          3          4
 1      2.0396    -1.5297     0.2746    -0.0490
 2                 1.6401    -0.2500     0.6737
 3                            0.7887     0.6617
 4                                       0.5347

PDF version (NAG web site, 64-bit version, 64-bit version)
Chapter Contents
Chapter Introduction
NAG Toolbox

© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015