doublePOsolve(3) | LAPACK | doublePOsolve(3) |
doublePOsolve - double
subroutine dposv (UPLO, N, NRHS, A, LDA, B, LDB, INFO)
DPOSV computes the solution to system of linear equations A * X = B for PO
matrices subroutine dposvx (FACT, UPLO, N, NRHS, A, LDA, AF,
LDAF, EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO)
DPOSVX computes the solution to system of linear equations A * X = B for
PO matrices subroutine dposvxx (FACT, UPLO, N, NRHS, A, LDA, AF,
LDAF, EQUED, S, B, LDB, X, LDX, RCOND, RPVGRW, BERR, N_ERR_BNDS,
ERR_BNDS_NORM, ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, IWORK, INFO)
DPOSVXX computes the solution to system of linear equations A * X = B for
PO matrices subroutine dsposv (UPLO, N, NRHS, A, LDA, B, LDB, X,
LDX, WORK, SWORK, ITER, INFO)
DSPOSV computes the solution to system of linear equations A * X = B for
PO matrices
This is the group of double solve driver functions for PO matrices
DPOSV computes the solution to system of linear equations A * X = B for PO matrices
Purpose:
DPOSV computes the solution to a real system of linear equations
A * X = B,
where A is an N-by-N symmetric positive definite matrix and X and B
are N-by-NRHS matrices.
The Cholesky decomposition is used to factor A as
A = U**T* U, if UPLO = 'U', or
A = L * L**T, 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 A * X = B.
Parameters
UPLO is CHARACTER*1
= 'U': Upper triangle of A is stored;
= 'L': Lower triangle of A is stored.
N
N is INTEGER
The number of linear equations, i.e., the order of the
matrix A. N >= 0.
NRHS
NRHS is INTEGER
The number of right hand sides, i.e., the number of columns
of the matrix B. NRHS >= 0.
A
A is DOUBLE PRECISION array, dimension (LDA,N)
On entry, the symmetric matrix A. If UPLO = 'U', the leading
N-by-N upper triangular part of A contains the upper
triangular part of the matrix A, and the strictly lower
triangular part of A is not referenced. If UPLO = 'L', the
leading N-by-N lower triangular part of A contains the lower
triangular part of the matrix A, and the strictly upper
triangular part of A is not referenced.
On exit, if INFO = 0, the factor U or L from the Cholesky
factorization A = U**T*U or A = L*L**T.
LDA
LDA is INTEGER
The leading dimension of the array A. LDA >= max(1,N).
B
B is DOUBLE PRECISION array, dimension (LDB,NRHS)
On entry, the N-by-NRHS right hand side matrix B.
On exit, if INFO = 0, the N-by-NRHS solution matrix X.
LDB
LDB is INTEGER
The leading dimension of the array B. LDB >= max(1,N).
INFO
INFO is INTEGER
= 0: successful exit
< 0: if INFO = -i, the i-th argument had an illegal value
> 0: if INFO = i, the leading minor of order i of A is not
positive definite, so the factorization could not be
completed, and the solution has not been computed.
Author
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
DPOSVX computes the solution to system of linear equations A * X = B for PO matrices
Purpose:
DPOSVX uses the Cholesky factorization A = U**T*U or A = L*L**T to
compute the solution to a real system of linear equations
A * X = B,
where A is an N-by-N symmetric positive definite matrix and X and B
are N-by-NRHS matrices.
Error bounds on the solution and a condition estimate are also
provided.
Description:
The following steps are performed:
1. If FACT = 'E', real scaling factors are computed to equilibrate
the system:
diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * B
Whether or not the system will be equilibrated depends on the
scaling of the matrix A, but if equilibration is used, A is
overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
factor the matrix A (after equilibration if FACT = 'E') as
A = U**T* U, if UPLO = 'U', or
A = L * L**T, if UPLO = 'L',
where U is an upper triangular matrix and L is a lower triangular
matrix.
3. If the leading i-by-i principal minor is not positive definite,
then the routine returns with INFO = i. Otherwise, the factored
form of A is used to estimate the condition number of the matrix
A. If the reciprocal of the condition number is less than machine
precision, INFO = N+1 is returned as a warning, but the routine
still goes on to solve for X and compute error bounds as
described below.
4. The system of equations is solved for X using the factored form
of A.
5. Iterative refinement is applied to improve the computed solution
matrix and calculate error bounds and backward error estimates
for it.
6. If equilibration was used, the matrix X is premultiplied by
diag(S) so that it solves the original system before
equilibration.
Parameters
FACT is CHARACTER*1
Specifies whether or not the factored form of the matrix A is
supplied on entry, and if not, whether the matrix A should be
equilibrated before it is factored.
= 'F': On entry, AF contains the factored form of A.
If EQUED = 'Y', the matrix A has been equilibrated
with scaling factors given by S. A and AF will not
be modified.
= 'N': The matrix A will be copied to AF and factored.
= 'E': The matrix A will be equilibrated if necessary, then
copied to AF and factored.
UPLO
UPLO is CHARACTER*1
= 'U': Upper triangle of A is stored;
= 'L': Lower triangle of A is stored.
N
N is INTEGER
The number of linear equations, i.e., the order of the
matrix A. N >= 0.
NRHS
NRHS is INTEGER
The number of right hand sides, i.e., the number of columns
of the matrices B and X. NRHS >= 0.
A
A is DOUBLE PRECISION array, dimension (LDA,N)
On entry, the symmetric matrix A, except if FACT = 'F' and
EQUED = 'Y', then A must contain the equilibrated matrix
diag(S)*A*diag(S). If UPLO = 'U', the leading
N-by-N upper triangular part of A contains the upper
triangular part of the matrix A, and the strictly lower
triangular part of A is not referenced. If UPLO = 'L', the
leading N-by-N lower triangular part of A contains the lower
triangular part of the matrix A, and the strictly upper
triangular part of A is not referenced. A is not modified if
FACT = 'F' or 'N', or if FACT = 'E' and EQUED = 'N' on exit.
On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
diag(S)*A*diag(S).
LDA
LDA is INTEGER
The leading dimension of the array A. LDA >= max(1,N).
AF
AF is DOUBLE PRECISION array, dimension (LDAF,N)
If FACT = 'F', then AF is an input argument and on entry
contains the triangular factor U or L from the Cholesky
factorization A = U**T*U or A = L*L**T, in the same storage
format as A. If EQUED .ne. 'N', then AF is the factored form
of the equilibrated matrix diag(S)*A*diag(S).
If FACT = 'N', then AF is an output argument and on exit
returns the triangular factor U or L from the Cholesky
factorization A = U**T*U or A = L*L**T of the original
matrix A.
If FACT = 'E', then AF is an output argument and on exit
returns the triangular factor U or L from the Cholesky
factorization A = U**T*U or A = L*L**T of the equilibrated
matrix A (see the description of A for the form of the
equilibrated matrix).
LDAF
LDAF is INTEGER
The leading dimension of the array AF. LDAF >= max(1,N).
EQUED
EQUED is CHARACTER*1
Specifies the form of equilibration that was done.
= 'N': No equilibration (always true if FACT = 'N').
= 'Y': Equilibration was done, i.e., A has been replaced by
diag(S) * A * diag(S).
EQUED is an input argument if FACT = 'F'; otherwise, it is an
output argument.
S
S is DOUBLE PRECISION array, dimension (N)
The scale factors for A; not accessed if EQUED = 'N'. S is
an input argument if FACT = 'F'; otherwise, S is an output
argument. If FACT = 'F' and EQUED = 'Y', each element of S
must be positive.
B
B is DOUBLE PRECISION array, dimension (LDB,NRHS)
On entry, the N-by-NRHS right hand side matrix B.
On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y',
B is overwritten by diag(S) * B.
LDB
LDB is INTEGER
The leading dimension of the array B. LDB >= max(1,N).
X
X is DOUBLE PRECISION array, dimension (LDX,NRHS)
If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X to
the original system of equations. Note that if EQUED = 'Y',
A and B are modified on exit, and the solution to the
equilibrated system is inv(diag(S))*X.
LDX
LDX is INTEGER
The leading dimension of the array X. LDX >= max(1,N).
RCOND
RCOND is DOUBLE PRECISION
The estimate of the reciprocal condition number of the matrix
A after equilibration (if done). If RCOND is less than the
machine precision (in particular, if RCOND = 0), the matrix
is singular to working precision. This condition is
indicated by a return code of INFO > 0.
FERR
FERR is DOUBLE PRECISION array, dimension (NRHS)
The estimated forward error bound for each solution vector
X(j) (the j-th column of the solution matrix X).
If XTRUE is the true solution corresponding to X(j), FERR(j)
is an estimated upper bound for the magnitude of the largest
element in (X(j) - XTRUE) divided by the magnitude of the
largest element in X(j). The estimate is as reliable as
the estimate for RCOND, and is almost always a slight
overestimate of the true error.
BERR
BERR is DOUBLE PRECISION array, dimension (NRHS)
The componentwise relative backward error of each solution
vector X(j) (i.e., the smallest relative change in
any element of A or B that makes X(j) an exact solution).
WORK
WORK is DOUBLE PRECISION array, dimension (3*N)
IWORK
IWORK is INTEGER array, dimension (N)
INFO
INFO is INTEGER
= 0: successful exit
< 0: if INFO = -i, the i-th argument had an illegal value
> 0: if INFO = i, and i is
<= N: the leading minor of order i of A is
not positive definite, so the factorization
could not be completed, and the solution has not
been computed. RCOND = 0 is returned.
= N+1: U 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.
Author
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
DPOSVXX computes the solution to system of linear equations A * X = B for PO matrices
Purpose:
DPOSVXX uses the Cholesky factorization A = U**T*U or A = L*L**T
to compute the solution to a double precision system of linear equations
A * X = B, where A is an N-by-N symmetric positive definite matrix
and X and B are N-by-NRHS matrices.
If requested, both normwise and maximum componentwise error bounds
are returned. DPOSVXX will return a solution with a tiny
guaranteed error (O(eps) where eps is the working machine
precision) unless the matrix is very ill-conditioned, in which
case a warning is returned. Relevant condition numbers also are
calculated and returned.
DPOSVXX accepts user-provided factorizations and equilibration
factors; see the definitions of the FACT and EQUED options.
Solving with refinement and using a factorization from a previous
DPOSVXX call will also produce a solution with either O(eps)
errors or warnings, but we cannot make that claim for general
user-provided factorizations and equilibration factors if they
differ from what DPOSVXX would itself produce.
Description:
The following steps are performed:
1. If FACT = 'E', double precision scaling factors are computed to equilibrate
the system:
diag(S)*A*diag(S) *inv(diag(S))*X = diag(S)*B
Whether or not the system will be equilibrated depends on the
scaling of the matrix A, but if equilibration is used, A is
overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
factor the matrix A (after equilibration if FACT = 'E') as
A = U**T* U, if UPLO = 'U', or
A = L * L**T, if UPLO = 'L',
where U is an upper triangular matrix and L is a lower triangular
matrix.
3. If the leading i-by-i principal minor is not positive definite,
then the routine returns with INFO = i. Otherwise, the factored
form of A is used to estimate the condition number of the matrix
A (see argument RCOND). If the reciprocal of the condition number
is less than machine precision, the routine still goes on to solve
for X and compute error bounds as described below.
4. The system of equations is solved for X using the factored form
of A.
5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
the routine will use iterative refinement to try to get a small
error and error bounds. Refinement calculates the residual to at
least twice the working precision.
6. If equilibration was used, the matrix X is premultiplied by
diag(S) so that it solves the original system before
equilibration.
Some optional parameters are bundled in the PARAMS array. These
settings determine how refinement is performed, but often the
defaults are acceptable. If the defaults are acceptable, users
can pass NPARAMS = 0 which prevents the source code from accessing
the PARAMS argument.
Parameters
FACT is CHARACTER*1
Specifies whether or not the factored form of the matrix A is
supplied on entry, and if not, whether the matrix A should be
equilibrated before it is factored.
= 'F': On entry, AF contains the factored form of A.
If EQUED is not 'N', the matrix A has been
equilibrated with scaling factors given by S.
A and AF are not modified.
= 'N': The matrix A will be copied to AF and factored.
= 'E': The matrix A will be equilibrated if necessary, then
copied to AF and factored.
UPLO
UPLO is CHARACTER*1
= 'U': Upper triangle of A is stored;
= 'L': Lower triangle of A is stored.
N
N is INTEGER
The number of linear equations, i.e., the order of the
matrix A. N >= 0.
NRHS
NRHS is INTEGER
The number of right hand sides, i.e., the number of columns
of the matrices B and X. NRHS >= 0.
A
A is DOUBLE PRECISION array, dimension (LDA,N)
On entry, the symmetric matrix A, except if FACT = 'F' and EQUED =
'Y', then A must contain the equilibrated matrix
diag(S)*A*diag(S). If UPLO = 'U', the leading N-by-N upper
triangular part of A contains the upper triangular part of the
matrix A, and the strictly lower triangular part of A is not
referenced. If UPLO = 'L', the leading N-by-N lower triangular
part of A contains the lower triangular part of the matrix A, and
the strictly upper triangular part of A is not referenced. A is
not modified if FACT = 'F' or 'N', or if FACT = 'E' and EQUED =
'N' on exit.
On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
diag(S)*A*diag(S).
LDA
LDA is INTEGER
The leading dimension of the array A. LDA >= max(1,N).
AF
AF is DOUBLE PRECISION array, dimension (LDAF,N)
If FACT = 'F', then AF is an input argument and on entry
contains the triangular factor U or L from the Cholesky
factorization A = U**T*U or A = L*L**T, in the same storage
format as A. If EQUED .ne. 'N', then AF is the factored
form of the equilibrated matrix diag(S)*A*diag(S).
If FACT = 'N', then AF is an output argument and on exit
returns the triangular factor U or L from the Cholesky
factorization A = U**T*U or A = L*L**T of the original
matrix A.
If FACT = 'E', then AF is an output argument and on exit
returns the triangular factor U or L from the Cholesky
factorization A = U**T*U or A = L*L**T of the equilibrated
matrix A (see the description of A for the form of the
equilibrated matrix).
LDAF
LDAF is INTEGER
The leading dimension of the array AF. LDAF >= max(1,N).
EQUED
EQUED is CHARACTER*1
Specifies the form of equilibration that was done.
= 'N': No equilibration (always true if FACT = 'N').
= 'Y': Both row and column equilibration, i.e., A has been
replaced by diag(S) * A * diag(S).
EQUED is an input argument if FACT = 'F'; otherwise, it is an
output argument.
S
S is DOUBLE PRECISION array, dimension (N)
The row scale factors for A. If EQUED = 'Y', A is multiplied on
the left and right by diag(S). S is an input argument if FACT =
'F'; otherwise, S is an output argument. If FACT = 'F' and EQUED
= 'Y', each element of S must be positive. If S is output, each
element of S is a power of the radix. If S is input, each element
of S should be a power of the radix to ensure a reliable solution
and error estimates. Scaling by powers of the radix does not cause
rounding errors unless the result underflows or overflows.
Rounding errors during scaling lead to refining with a matrix that
is not equivalent to the input matrix, producing error estimates
that may not be reliable.
B
B is DOUBLE PRECISION array, dimension (LDB,NRHS)
On entry, the N-by-NRHS right hand side matrix B.
On exit,
if EQUED = 'N', B is not modified;
if EQUED = 'Y', B is overwritten by diag(S)*B;
LDB
LDB is INTEGER
The leading dimension of the array B. LDB >= max(1,N).
X
X is DOUBLE PRECISION array, dimension (LDX,NRHS)
If INFO = 0, the N-by-NRHS solution matrix X to the original
system of equations. Note that A and B are modified on exit if
EQUED .ne. 'N', and the solution to the equilibrated system is
inv(diag(S))*X.
LDX
LDX is INTEGER
The leading dimension of the array X. LDX >= max(1,N).
RCOND
RCOND is DOUBLE PRECISION
Reciprocal scaled condition number. This is an estimate of the
reciprocal Skeel condition number of the matrix A after
equilibration (if done). If this is less than the machine
precision (in particular, if it is zero), the matrix is singular
to working precision. Note that the error may still be small even
if this number is very small and the matrix appears ill-
conditioned.
RPVGRW
RPVGRW is DOUBLE PRECISION
Reciprocal pivot growth. On exit, this contains the reciprocal
pivot growth factor norm(A)/norm(U). The 'max absolute element'
norm is used. If this is much less than 1, then the stability of
the LU factorization of the (equilibrated) matrix A could be poor.
This also means that the solution X, estimated condition numbers,
and error bounds could be unreliable. If factorization fails with
0<INFO<=N, then this contains the reciprocal pivot growth factor
for the leading INFO columns of A.
BERR
BERR is DOUBLE PRECISION array, dimension (NRHS)
Componentwise relative backward error. This is the
componentwise relative backward error of each solution vector X(j)
(i.e., the smallest relative change in any element of A or B that
makes X(j) an exact solution).
N_ERR_BNDS
N_ERR_BNDS is INTEGER
Number of error bounds to return for each right hand side
and each type (normwise or componentwise). See ERR_BNDS_NORM and
ERR_BNDS_COMP below.
ERR_BNDS_NORM
ERR_BNDS_NORM is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
For each right-hand side, this array contains information about
various error bounds and condition numbers corresponding to the
normwise relative error, which is defined as follows:
Normwise relative error in the ith solution vector:
max_j (abs(XTRUE(j,i) - X(j,i)))
------------------------------
max_j abs(X(j,i))
The array is indexed by the type of error information as described
below. There currently are up to three pieces of information
returned.
The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
right-hand side.
The second index in ERR_BNDS_NORM(:,err) contains the following
three fields:
err = 1 'Trust/don't trust' boolean. Trust the answer if the
reciprocal condition number is less than the threshold
sqrt(n) * dlamch('Epsilon').
err = 2 'Guaranteed' error bound: The estimated forward error,
almost certainly within a factor of 10 of the true error
so long as the next entry is greater than the threshold
sqrt(n) * dlamch('Epsilon'). This error bound should only
be trusted if the previous boolean is true.
err = 3 Reciprocal condition number: Estimated normwise
reciprocal condition number. Compared with the threshold
sqrt(n) * dlamch('Epsilon') to determine if the error
estimate is 'guaranteed'. These reciprocal condition
numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
appropriately scaled matrix Z.
Let Z = S*A, where S scales each row by a power of the
radix so all absolute row sums of Z are approximately 1.
See Lapack Working Note 165 for further details and extra
cautions.
ERR_BNDS_COMP
ERR_BNDS_COMP is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
For each right-hand side, this array contains information about
various error bounds and condition numbers corresponding to the
componentwise relative error, which is defined as follows:
Componentwise relative error in the ith solution vector:
abs(XTRUE(j,i) - X(j,i))
max_j ----------------------
abs(X(j,i))
The array is indexed by the right-hand side i (on which the
componentwise relative error depends), and the type of error
information as described below. There currently are up to three
pieces of information returned for each right-hand side. If
componentwise accuracy is not requested (PARAMS(3) = 0.0), then
ERR_BNDS_COMP is not accessed. If N_ERR_BNDS < 3, then at most
the first (:,N_ERR_BNDS) entries are returned.
The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
right-hand side.
The second index in ERR_BNDS_COMP(:,err) contains the following
three fields:
err = 1 'Trust/don't trust' boolean. Trust the answer if the
reciprocal condition number is less than the threshold
sqrt(n) * dlamch('Epsilon').
err = 2 'Guaranteed' error bound: The estimated forward error,
almost certainly within a factor of 10 of the true error
so long as the next entry is greater than the threshold
sqrt(n) * dlamch('Epsilon'). This error bound should only
be trusted if the previous boolean is true.
err = 3 Reciprocal condition number: Estimated componentwise
reciprocal condition number. Compared with the threshold
sqrt(n) * dlamch('Epsilon') to determine if the error
estimate is 'guaranteed'. These reciprocal condition
numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
appropriately scaled matrix Z.
Let Z = S*(A*diag(x)), where x is the solution for the
current right-hand side and S scales each row of
A*diag(x) by a power of the radix so all absolute row
sums of Z are approximately 1.
See Lapack Working Note 165 for further details and extra
cautions.
NPARAMS
NPARAMS is INTEGER
Specifies the number of parameters set in PARAMS. If <= 0, the
PARAMS array is never referenced and default values are used.
PARAMS
PARAMS is DOUBLE PRECISION array, dimension NPARAMS
Specifies algorithm parameters. If an entry is < 0.0, then
that entry will be filled with default value used for that
parameter. Only positions up to NPARAMS are accessed; defaults
are used for higher-numbered parameters.
PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
refinement or not.
Default: 1.0D+0
= 0.0: No refinement is performed, and no error bounds are
computed.
= 1.0: Use the extra-precise refinement algorithm.
(other values are reserved for future use)
PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
computations allowed for refinement.
Default: 10
Aggressive: Set to 100 to permit convergence using approximate
factorizations or factorizations other than LU. If
the factorization uses a technique other than
Gaussian elimination, the guarantees in
err_bnds_norm and err_bnds_comp may no longer be
trustworthy.
PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
will attempt to find a solution with small componentwise
relative error in the double-precision algorithm. Positive
is true, 0.0 is false.
Default: 1.0 (attempt componentwise convergence)
WORK
WORK is DOUBLE PRECISION array, dimension (4*N)
IWORK
IWORK is INTEGER array, dimension (N)
INFO
INFO is INTEGER
= 0: Successful exit. The solution to every right-hand side is
guaranteed.
< 0: If INFO = -i, the i-th argument had an illegal value
> 0 and <= N: U(INFO,INFO) is exactly zero. The factorization
has been completed, but the factor U is exactly singular, so
the solution and error bounds could not be computed. RCOND = 0
is returned.
= N+J: The solution corresponding to the Jth right-hand side is
not guaranteed. The solutions corresponding to other right-
hand sides K with K > J may not be guaranteed as well, but
only the first such right-hand side is reported. If a small
componentwise error is not requested (PARAMS(3) = 0.0) then
the Jth right-hand side is the first with a normwise error
bound that is not guaranteed (the smallest J such
that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
the Jth right-hand side is the first with either a normwise or
componentwise error bound that is not guaranteed (the smallest
J such that either ERR_BNDS_NORM(J,1) = 0.0 or
ERR_BNDS_COMP(J,1) = 0.0). See the definition of
ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
about all of the right-hand sides check ERR_BNDS_NORM or
ERR_BNDS_COMP.
Author
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
DSPOSV computes the solution to system of linear equations A * X = B for PO matrices
Purpose:
DSPOSV computes the solution to a real system of linear equations
A * X = B,
where A is an N-by-N symmetric positive definite matrix and X and B
are N-by-NRHS matrices.
DSPOSV first attempts to factorize the matrix in SINGLE PRECISION
and use this factorization within an iterative refinement procedure
to produce a solution with DOUBLE PRECISION normwise backward error
quality (see below). If the approach fails the method switches to a
DOUBLE PRECISION factorization and solve.
The iterative refinement is not going to be a winning strategy if
the ratio SINGLE PRECISION performance over DOUBLE PRECISION
performance is too small. A reasonable strategy should take the
number of right-hand sides and the size of the matrix into account.
This might be done with a call to ILAENV in the future. Up to now, we
always try iterative refinement.
The iterative refinement process is stopped if
ITER > ITERMAX
or for all the RHS we have:
RNRM < SQRT(N)*XNRM*ANRM*EPS*BWDMAX
where
o ITER is the number of the current iteration in the iterative
refinement process
o RNRM is the infinity-norm of the residual
o XNRM is the infinity-norm of the solution
o ANRM is the infinity-operator-norm of the matrix A
o EPS is the machine epsilon returned by DLAMCH('Epsilon')
The value ITERMAX and BWDMAX are fixed to 30 and 1.0D+00
respectively.
Parameters
UPLO is CHARACTER*1
= 'U': Upper triangle of A is stored;
= 'L': Lower triangle of A is stored.
N
N is INTEGER
The number of linear equations, i.e., the order of the
matrix A. N >= 0.
NRHS
NRHS is INTEGER
The number of right hand sides, i.e., the number of columns
of the matrix B. NRHS >= 0.
A
A is DOUBLE PRECISION array,
dimension (LDA,N)
On entry, the symmetric matrix A. If UPLO = 'U', the leading
N-by-N upper triangular part of A contains the upper
triangular part of the matrix A, and the strictly lower
triangular part of A is not referenced. If UPLO = 'L', the
leading N-by-N lower triangular part of A contains the lower
triangular part of the matrix A, and the strictly upper
triangular part of A is not referenced.
On exit, if iterative refinement has been successfully used
(INFO = 0 and ITER >= 0, see description below), then A is
unchanged, if double precision factorization has been used
(INFO = 0 and ITER < 0, see description below), then the
array A contains the factor U or L from the Cholesky
factorization A = U**T*U or A = L*L**T.
LDA
LDA is INTEGER
The leading dimension of the array A. LDA >= max(1,N).
B
B is DOUBLE PRECISION array, dimension (LDB,NRHS)
The N-by-NRHS right hand side matrix B.
LDB
LDB is INTEGER
The leading dimension of the array B. LDB >= max(1,N).
X
X is DOUBLE PRECISION array, dimension (LDX,NRHS)
If INFO = 0, the N-by-NRHS solution matrix X.
LDX
LDX is INTEGER
The leading dimension of the array X. LDX >= max(1,N).
WORK
WORK is DOUBLE PRECISION array, dimension (N,NRHS)
This array is used to hold the residual vectors.
SWORK
SWORK is REAL array, dimension (N*(N+NRHS))
This array is used to use the single precision matrix and the
right-hand sides or solutions in single precision.
ITER
ITER is INTEGER
< 0: iterative refinement has failed, double precision
factorization has been performed
-1 : the routine fell back to full precision for
implementation- or machine-specific reasons
-2 : narrowing the precision induced an overflow,
the routine fell back to full precision
-3 : failure of SPOTRF
-31: stop the iterative refinement after the 30th
iterations
> 0: iterative refinement has been successfully used.
Returns the number of iterations
INFO
INFO is INTEGER
= 0: successful exit
< 0: if INFO = -i, the i-th argument had an illegal value
> 0: if INFO = i, the leading minor of order i of (DOUBLE
PRECISION) A is not positive definite, so the
factorization could not be completed, and the solution
has not been computed.
Author
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Generated automatically by Doxygen for LAPACK from the source code.
Sun Nov 27 2022 | Version 3.11.0 |