NAG FL Interface
f07hbf (dpbsvx)
1
Purpose
f07hbf uses the Cholesky factorization
to compute the solution to a real system of linear equations
where
is an
by
symmetric positive definite band matrix of bandwidth
and
and
are
by
matrices. Error bounds on the solution and a condition estimate are also provided.
2
Specification
Fortran Interface
Subroutine f07hbf ( |
fact, uplo, n, kd, nrhs, ab, ldab, afb, ldafb, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, iwork, info) |
Integer, Intent (In) |
:: |
n, kd, nrhs, ldab, ldafb, ldb, ldx |
Integer, Intent (Out) |
:: |
iwork(n), info |
Real (Kind=nag_wp), Intent (Inout) |
:: |
ab(ldab,*), afb(ldafb,*), s(*), b(ldb,*), x(ldx,*) |
Real (Kind=nag_wp), Intent (Out) |
:: |
rcond, ferr(nrhs), berr(nrhs), work(3*n) |
Character (1), Intent (In) |
:: |
fact, uplo |
Character (1), Intent (InOut) |
:: |
equed |
|
C Header Interface
#include <nag.h>
void |
f07hbf_ (const char *fact, const char *uplo, const Integer *n, const Integer *kd, const Integer *nrhs, double ab[], const Integer *ldab, double afb[], const Integer *ldafb, char *equed, double s[], double b[], const Integer *ldb, double x[], const Integer *ldx, double *rcond, double ferr[], double berr[], double work[], Integer iwork[], Integer *info, const Charlen length_fact, const Charlen length_uplo, const Charlen length_equed) |
|
C++ Header Interface
#include <nag.h> extern "C" {
void |
f07hbf_ (const char *fact, const char *uplo, const Integer &n, const Integer &kd, const Integer &nrhs, double ab[], const Integer &ldab, double afb[], const Integer &ldafb, char *equed, double s[], double b[], const Integer &ldb, double x[], const Integer &ldx, double &rcond, double ferr[], double berr[], double work[], Integer iwork[], Integer &info, const Charlen length_fact, const Charlen length_uplo, const Charlen length_equed) |
}
|
The routine may be called by the names f07hbf, nagf_lapacklin_dpbsvx or its LAPACK name dpbsvx.
3
Description
f07hbf performs the following steps:
-
1.If , real diagonal scaling factors, , are computed to equilibrate the system:
Whether or not the system will be equilibrated depends on the scaling of the matrix , but if equilibration is used, is overwritten by and by .
-
2.If or , the Cholesky decomposition is used to factor the matrix (after equilibration if ) as if or if , where is an upper triangular matrix and is a lower triangular matrix.
-
3.If the leading by principal minor of is not positive definite, then the routine returns with . Otherwise, the factored form of is used to estimate the condition number of the matrix . If the reciprocal of the condition number is less than machine precision, is returned as a warning, but the routine still goes on to solve for and compute error bounds as described below.
-
4.The system of equations is solved for using the factored form of .
-
5.Iterative refinement is applied to improve the computed solution matrix and to calculate error bounds and backward error estimates for it.
-
6.If equilibration was used, the matrix is premultiplied by so that it solves the original system before equilibration.
4
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
https://www.netlib.org/lapack/lug
Golub G H and Van Loan C F (1996) Matrix Computations (3rd Edition) Johns Hopkins University Press, Baltimore
Higham N J (2002) Accuracy and Stability of Numerical Algorithms (2nd Edition) SIAM, Philadelphia
5
Arguments
-
1:
– Character(1)
Input
-
On entry: specifies whether or not the factorized form of the matrix
is supplied on entry, and if not, whether the matrix
should be equilibrated before it is factorized.
- afb contains the factorized form of . If , the matrix has been equilibrated with scaling factors given by s. ab and afb will not be modified.
- The matrix will be copied to afb and factorized.
- The matrix will be equilibrated if necessary, then copied to afb and factorized.
Constraint:
, or .
-
2:
– Character(1)
Input
-
On entry: if
, the upper triangle of
is stored.
If , the lower triangle of is stored.
Constraint:
or .
-
3:
– Integer
Input
-
On entry: , the number of linear equations, i.e., the order of the matrix .
Constraint:
.
-
4:
– Integer
Input
-
On entry: , the number of superdiagonals of the matrix if , or the number of subdiagonals if .
Constraint:
.
-
5:
– Integer
Input
-
On entry: , the number of right-hand sides, i.e., the number of columns of the matrix .
Constraint:
.
-
6:
– Real (Kind=nag_wp) array
Input/Output
-
Note: the second dimension of the array
ab
must be at least
.
On entry: the upper or lower triangle of the symmetric band matrix
, except if
and
, in which case
ab must contain the equilibrated matrix
.
The matrix is stored in rows
to
, more precisely,
- if , the elements of the upper triangle of within the band must be stored with element in ;
- if , the elements of the lower triangle of within the band must be stored with element in
On exit: if
and
,
ab is overwritten by
.
-
7:
– Integer
Input
-
On entry: the first dimension of the array
ab as declared in the (sub)program from which
f07hbf is called.
Constraint:
.
-
8:
– Real (Kind=nag_wp) array
Input/Output
-
Note: the second dimension of the array
afb
must be at least
.
On entry: if
,
afb contains the triangular factor
or
from the Cholesky factorization
or
of the band matrix
, in the same storage format as
. If
,
afb is the factorized form of the equilibrated matrix
.
On exit: if
,
afb returns the triangular factor
or
from the Cholesky factorization
or
.
If
,
afb returns the triangular factor
or
from the Cholesky factorization
or
of the equilibrated matrix
(see the description of
ab for the form of the equilibrated matrix).
-
9:
– Integer
Input
-
On entry: the first dimension of the array
afb as declared in the (sub)program from which
f07hbf is called.
Constraint:
.
-
10:
– Character(1)
Input/Output
-
On entry: if
or
,
equed need not be set.
If
,
equed must specify the form of the equilibration that was performed as follows:
- if , no equilibration;
- if , equilibration was performed, i.e., has been replaced by .
On exit: if
,
equed is unchanged from entry.
Otherwise, if no constraints are violated,
equed specifies the form of the equilibration that was performed as specified above.
Constraint:
if , or .
-
11:
– Real (Kind=nag_wp) array
Input/Output
-
Note: the dimension of the array
s
must be at least
.
On entry: if
or
,
s need not be set.
If
and
,
s must contain the scale factors,
, for
; each element of
s must be positive.
On exit: if
,
s is unchanged from entry.
Otherwise, if no constraints are violated and
,
s contains the scale factors,
, for
; each element of
s is positive.
-
12:
– Real (Kind=nag_wp) array
Input/Output
-
Note: the second dimension of the array
b
must be at least
.
On entry: the by right-hand side matrix .
On exit: if
,
b is not modified.
If
,
b is overwritten by
.
-
13:
– Integer
Input
-
On entry: the first dimension of the array
b as declared in the (sub)program from which
f07hbf is called.
Constraint:
.
-
14:
– Real (Kind=nag_wp) array
Output
-
Note: the second dimension of the array
x
must be at least
.
On exit: if or , the by solution matrix to the original system of equations. Note that the arrays and are modified on exit if , and the solution to the equilibrated system is .
-
15:
– Integer
Input
-
On entry: the first dimension of the array
x as declared in the (sub)program from which
f07hbf is called.
Constraint:
.
-
16:
– Real (Kind=nag_wp)
Output
-
On exit: if no constraints are violated, an estimate of the reciprocal condition number of the matrix (after equilibration if that is performed), computed as .
-
17:
– Real (Kind=nag_wp) array
Output
-
On exit: if
or
, an estimate of the forward error bound for each computed solution vector, such that
where
is the
th column of the computed solution returned in the array
x and
is the corresponding column of the exact solution
. The estimate is as reliable as the estimate for
rcond, and is almost always a slight overestimate of the true error.
-
18:
– Real (Kind=nag_wp) array
Output
-
On exit: if or , an estimate of the component-wise relative backward error of each computed solution vector (i.e., the smallest relative change in any element of or that makes an exact solution).
-
19:
– Real (Kind=nag_wp) array
Workspace
-
-
20:
– Integer array
Workspace
-
-
21:
– Integer
Output
-
On exit:
unless the routine detects an error (see
Section 6).
6
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. is returned.
-
(or
) is nonsingular, but
rcond is less than
machine precision, meaning that the matrix is singular to working precision.
Nevertheless, the solution and error bounds are computed because there
are a number of situations where the computed solution can be more accurate
than the value of
rcond would suggest.
7
Accuracy
For each right-hand side vector
, the computed solution
is the exact solution of a perturbed system of equations
, where
- if , ;
- if , ,
is a modest linear function of
, and
is the
machine precision. See Section 10.1 of
Higham (2002) for further details.
If
is the true solution, then the computed solution
satisfies a forward error bound of the form
where
.
If
is the
th column of
, then
is returned in
and a bound on
is returned in
. See Section 4.4 of
Anderson et al. (1999) for further details.
8
Parallelism and Performance
f07hbf is threaded by NAG for parallel execution in multithreaded implementations of the NAG Library.
f07hbf makes calls to BLAS and/or LAPACK routines, which may be threaded within the vendor library used by this implementation. Consult the documentation for the vendor library for further information.
Please consult the
X06 Chapter Introduction for information on how to control and interrogate the OpenMP environment used within this routine. Please also consult the
Users' Note for your implementation for any additional implementation-specific information.
When , the factorization of requires approximately floating-point operations, where is the number of superdiagonals.
For each right-hand side, computation of the backward error involves a minimum of floating-point operations. Each step of iterative refinement involves an additional operations. At most five steps of iterative refinement are performed, but usually only one or two steps are required. Estimating the forward error involves solving a number of systems of equations of the form ; the number is usually or and never more than . Each solution involves approximately operations.
The complex analogue of this routine is
f07hpf.
10
Example
This example solves the equations
where
is the symmetric positive definite band matrix
and
Error estimates for the solutions, information on equilibration and an estimate of the reciprocal of the condition number of the scaled matrix are also output.
10.1
Program Text
10.2
Program Data
10.3
Program Results