PDF version (NAG web site
, 64-bit version, 64-bit version)
NAG Toolbox: nag_lapack_zppsv (f07gn)
Purpose
nag_lapack_zppsv (f07gn) computes the solution to a complex system of linear equations
where
is an
by
Hermitian positive definite matrix stored in packed format and
and
are
by
matrices.
Syntax
Description
nag_lapack_zppsv (f07gn) uses the Cholesky decomposition to factor as if or if , where is an upper triangular matrix and is a lower triangular matrix. The factored form of is then used to solve the system of equations .
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:
– string (length ≥ 1)
-
If
, the upper triangle of
is stored.
If , the lower triangle of is stored.
Constraint:
or .
- 2:
– complex array
-
The dimension of the array
ap
must be at least
The
by
Hermitian matrix
, packed by columns.
More precisely,
- if , the upper triangle of must be stored with element in for ;
- if , the lower triangle of must be stored with element in for .
- 3:
– complex array
-
The first dimension of the array
b must be at least
.
The second dimension of the array
b must be at least
.
Note: to solve the equations
, where
is a single right-hand side,
b may be supplied as a one-dimensional array with length
.
The by right-hand side matrix .
Optional Input Parameters
- 1:
– int64int32nag_int scalar
-
Default:
the first dimension of the array
b.
, the number of linear equations, i.e., the order of the matrix .
Constraint:
.
- 2:
– int64int32nag_int scalar
-
Default:
the second dimension of the array
b.
, the number of right-hand sides, i.e., the number of columns of the matrix .
Constraint:
.
Output Parameters
- 1:
– complex array
-
The dimension of the array
ap will be
If , the factor or from the Cholesky factorization or , in the same storage format as .
- 2:
– complex array
-
The first dimension of the array
b will be
.
The second dimension of the array
b will be
.
Note: to solve the equations
, where
is a single right-hand side,
b may be supplied as a one-dimensional array with length
.
If , the by solution matrix .
- 3:
– int64int32nag_int scalar
unless the function detects an error (see
Error Indicators and Warnings).
Error Indicators and Warnings
-
If , argument had an illegal value. An explanatory message is output, and execution of the program is terminated.
-
-
The leading minor of order of 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,
, satisfies an equation of the form
where
and
is the
machine precision. An approximate error bound for the computed solution is given by
where
, the condition number of
with respect to the solution of the linear equations. See Section 4.4 of
Anderson et al. (1999) for further details.
nag_lapack_zppsvx (f07gp) is a comprehensive LAPACK driver that returns forward and backward error bounds and an estimate of the condition number. Alternatively,
nag_linsys_complex_posdef_packed_solve (f04ce) solves
and returns a forward error bound and condition estimate.
nag_linsys_complex_posdef_packed_solve (f04ce) calls
nag_lapack_zppsv (f07gn) to solve the equations.
Further Comments
The total number of floating-point operations is approximately , where is the number of right-hand sides.
The real analogue of this function is
nag_lapack_dppsv (f07ga).
Example
This example solves the equations
where
is the Hermitian positive definite matrix
and
Details of the Cholesky factorization of are also output.
Open in the MATLAB editor:
f07gn_example
function f07gn_example
fprintf('f07gn example results\n\n');
uplo = 'Upper';
n = int64(4);
ap = [ 3.23 + 0i, ...
1.51 - 1.92i, 3.58 + 0i, ...
1.90 + 0.84i, -0.23 + 1.11i, 4.09 + 0i, ...
0.42 + 2.50i, -1.18 + 1.37i, 2.33 - 0.14i, 4.29 + 0i];
b = [ 3.93 - 6.14i;
6.17 + 9.42i;
-7.17 - 21.83i;
1.99 - 14.38i];
[apf, x, info] = f07gn( ...
uplo, ap, b);
disp('Solution');
disp(x);
[ifail] = x04dc( ...
uplo, 'Non-unit', n, apf, 'Cholesky factor U');
f07gn example results
Solution
1.0000 - 1.0000i
-0.0000 + 3.0000i
-4.0000 - 5.0000i
2.0000 + 1.0000i
Cholesky factor U
1 2 3 4
1 1.7972 0.8402 1.0572 0.2337
0.0000 -1.0683 0.4674 1.3910
2 1.3164 -0.4702 0.0834
0.0000 -0.3131 -0.0368
3 1.5604 0.9360
0.0000 -0.9900
4 0.6603
0.0000
PDF version (NAG web site
, 64-bit version, 64-bit version)
© The Numerical Algorithms Group Ltd, Oxford, UK. 2009–2015