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_zhetrs (f07ms)

 Contents

    1  Purpose
    2  Syntax
    7  Accuracy
    9  Example

Purpose

nag_lapack_zhetrs (f07ms) solves a complex Hermitian indefinite system of linear equations with multiple right-hand sides,
AX=B ,  
where A has been factorized by nag_lapack_zhetrf (f07mr).

Syntax

[b, info] = f07ms(uplo, a, ipiv, b, 'n', n, 'nrhs_p', nrhs_p)
[b, info] = nag_lapack_zhetrs(uplo, a, ipiv, b, 'n', n, 'nrhs_p', nrhs_p)

Description

nag_lapack_zhetrs (f07ms) is used to solve a complex Hermitian indefinite system of linear equations AX=B, this function must be preceded by a call to nag_lapack_zhetrf (f07mr) which computes the Bunch–Kaufman factorization of A.
If uplo='U', A=PUDUHPT, where P is a permutation matrix, U is an upper triangular matrix and D is an Hermitian block diagonal matrix with 1 by 1 and 2 by 2 blocks; the solution X is computed by solving PUDY=B and then UHPTX=Y.
If uplo='L', A=PLDLHPT, where L is a lower triangular matrix; the solution X is computed by solving PLDY=B and then LHPTX=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=PUDUHPT, where U is upper triangular.
uplo='L'
A=PLDLHPT, where L is lower triangular.
Constraint: uplo='U' or 'L'.
2:     alda: – complex array
The first dimension of the array a must be at least max1,n.
The second dimension of the array a must be at least max1,n.
Details of the factorization of A, as returned by nag_lapack_zhetrf (f07mr).
3:     ipiv: int64int32nag_int array
The dimension of the array ipiv must be at least max1,n
Details of the interchanges and the block structure of D, as returned by nag_lapack_zhetrf (f07mr).
4:     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 arrays a, b and the second dimension of the arrays a, ipiv.
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_zherfs (f07mv), and an estimate for κA (=κ1A) can be obtained by calling nag_lapack_zhecon (f07mu).

Further Comments

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

Example

This example solves the system of equations AX=B, where
A= -1.36+0.00i 1.58+0.90i 2.21-0.21i 3.91+1.50i 1.58-0.90i -8.87+0.00i -1.84-0.03i -1.78+1.18i 2.21+0.21i -1.84+0.03i -4.63+0.00i 0.11+0.11i 3.91-1.50i -1.78-1.18i 0.11-0.11i -1.84+0.00i  
and
B= 7.79+05.48i -35.39+18.01i -0.77-16.05i 4.23-70.02i -9.58+03.88i -24.79-08.40i 2.98-10.18i 28.68-39.89i .  
Here A is Hermitian indefinite and must first be factorized by nag_lapack_zhetrf (f07mr).
function f07ms_example


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

% Hermitian indefinite matrix A (Lower triangular part stored)
uplo = 'L';
a = [-1.36 + 0i,     0    + 0i,     0    + 0i,      0    + 0i;
      1.58 - 0.90i, -8.87 + 0i,     0    + 0i,      0    + 0i;
      2.21 + 0.21i, -1.84 + 0.03i, -4.63 + 0i,      0    + 0i;
      3.91 - 1.50i, -1.78 - 1.18i,  0.11 - 0.11i,  -1.84 + 0i];

% Factorize
[af, ipiv, info] = f07mr( ...
                          uplo, a);

% RHS
b = [ 7.79 +  5.48i, -35.39 + 18.01i;
     -0.77 - 16.05i,   4.23 - 70.02i;
     -9.58 +  3.88i, -24.79 -  8.40i;
      2.98 - 10.18i,  28.68 - 39.89i];

% Solve
[x, info] = f07ms( ...
                   uplo, af, ipiv, b);

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


f07ms example results

Solution(s)
   1.0000 - 1.0000i   3.0000 - 4.0000i
  -1.0000 + 2.0000i  -1.0000 + 5.0000i
   3.0000 - 2.0000i   7.0000 - 2.0000i
   2.0000 + 1.0000i  -8.0000 + 6.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