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_zpptrs (f07gs)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_lapack_zpptrs (f07gs) solves a complex Hermitian positive definite system of linear equations with multiple right-hand sides,
AX=B ,  
where A has been factorized by nag_lapack_zpptrf (f07gr), using packed storage.

Syntax

[b, info] = f07gs(uplo, ap, b, 'n', n, 'nrhs_p', nrhs_p)
[b, info] = nag_lapack_zpptrs(uplo, ap, b, 'n', n, 'nrhs_p', nrhs_p)

Description

nag_lapack_zpptrs (f07gs) is used to solve a complex Hermitian positive definite system of linear equations AX=B, the function must be preceded by a call to nag_lapack_zpptrf (f07gr) which computes the Cholesky factorization of A, using packed storage. The solution X is computed by forward and backward substitution.
If uplo='U', A=UHU, where U is upper triangular; the solution X is computed by solving UHY=B and then UX=Y.
If uplo='L', A=LLH, where L is lower triangular; the solution X is computed by solving LY=B and then LHX=Y.

References

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)
Specifies how A has been factorized.
uplo='U'
A=UHU, where U is upper triangular.
uplo='L'
A=LLH, where L is lower triangular.
Constraint: uplo='U' or 'L'.
2:     ap: – complex array
The dimension of the array ap must be at least max1,n×n+1/2
The Cholesky factor of A stored in packed form, as returned by nag_lapack_zpptrf (f07gr).
3:     bldb: – complex 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 ap and the second dimension of the array ap. (An error is raised if these dimensions are not equal.)
n, 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.
Constraint: nrhs_p0.

Output Parameters

1:     bldb: – complex array
The first dimension of the array b will be max1,n.
The second dimension of the array b will be max1,nrhs_p.
The n by r solution matrix X.
2:     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.

Accuracy

For each right-hand side vector b, the computed solution x is the exact solution of a perturbed system of equations A+Ex=b, where cn is a modest linear function of n, and ε is the machine precision
If x^ is the true solution, then the computed solution x satisfies a forward error bound of the form
x-x^ x cncondA,xε  
where condA,x=A-1Ax/xcondA=A-1AκA.
Note that condA,x can be much smaller than condA.
Forward and backward error bounds can be computed by calling nag_lapack_zpprfs (f07gv), and an estimate for κA (=κ1A) can be obtained by calling nag_lapack_zppcon (f07gu).

Further Comments

The total number of real floating-point operations is approximately 8n2r.
This function may be followed by a call to nag_lapack_zpprfs (f07gv) to refine the solution and return an error estimate.
The real analogue of this function is nag_lapack_dpptrs (f07ge).

Example

This example solves the system of equations AX=B, where
A= 3.23+0.00i 1.51-1.92i 1.90+0.84i 0.42+2.50i 1.51+1.92i 3.58+0.00i -0.23+1.11i -1.18+1.37i 1.90-0.84i -0.23-1.11i 4.09+0.00i 2.33-0.14i 0.42-2.50i -1.18-1.37i 2.33+0.14i 4.29+0.00i  
and
B= 3.93-06.14i 1.48+06.58i 6.17+09.42i 4.65-04.75i -7.17-21.83i -4.91+02.29i 1.99-14.38i 7.64-10.79i .  
Here A is Hermitian positive definite, stored in packed form, and must first be factorized by nag_lapack_zpptrf (f07gr).
function f07gs_example


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

% Symmetric A in lower packed form
uplo = 'L';
n = int64(4);
ap = [3.23 + 0i    1.51 + 1.92i    1.90 - 0.84i    0.42 - 2.50i ...
                   3.58 + 0i      -0.23 - 1.11i   -1.18 - 1.37i ...
                                   4.09 + 0.00i    2.33 + 0.14i ...
                                                   4.29 + 0.00i];

[L, info] = f07gr( ...
                     uplo, n, ap);

% RHS
b = [ 3.93 -  6.14i,  1.48 +  6.58i;
      6.17 +  9.42i,  4.65 -  4.75i;
     -7.17 - 21.83i, -4.91 +  2.29i;
      1.99 - 14.38i,  7.64 - 10.79i];

% Solve AX = B
  
[x, info] = f07gs( ...
                   uplo, L, b);

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


f07gs example results

Solution(s)
   1.0000 - 1.0000i  -1.0000 + 2.0000i
  -0.0000 + 3.0000i   3.0000 - 4.0000i
  -4.0000 - 5.0000i  -2.0000 + 3.0000i
   2.0000 + 1.0000i   4.0000 - 5.0000i


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