1 *> \brief <b> DPPSVX computes the solution to system of linear equations A * X = B for OTHER matrices</b>
3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download DPPSVX + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dppsvx.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dppsvx.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dppsvx.f">
21 * SUBROUTINE DPPSVX( FACT, UPLO, N, NRHS, AP, AFP, EQUED, S, B, LDB,
22 * X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO )
24 * .. Scalar Arguments ..
25 * CHARACTER EQUED, FACT, UPLO
26 * INTEGER INFO, LDB, LDX, N, NRHS
27 * DOUBLE PRECISION RCOND
29 * .. Array Arguments ..
31 * DOUBLE PRECISION AFP( * ), AP( * ), B( LDB, * ), BERR( * ),
32 * $ FERR( * ), S( * ), WORK( * ), X( LDX, * )
41 *> DPPSVX uses the Cholesky factorization A = U**T*U or A = L*L**T to
42 *> compute the solution to a real system of linear equations
44 *> where A is an N-by-N symmetric positive definite matrix stored in
45 *> packed format and X and B are N-by-NRHS matrices.
47 *> Error bounds on the solution and a condition estimate are also
56 *> The following steps are performed:
58 *> 1. If FACT = 'E', real scaling factors are computed to equilibrate
60 *> diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * B
61 *> Whether or not the system will be equilibrated depends on the
62 *> scaling of the matrix A, but if equilibration is used, A is
63 *> overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
65 *> 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
66 *> factor the matrix A (after equilibration if FACT = 'E') as
67 *> A = U**T* U, if UPLO = 'U', or
68 *> A = L * L**T, if UPLO = 'L',
69 *> where U is an upper triangular matrix and L is a lower triangular
72 *> 3. If the leading i-by-i principal minor is not positive definite,
73 *> then the routine returns with INFO = i. Otherwise, the factored
74 *> form of A is used to estimate the condition number of the matrix
75 *> A. If the reciprocal of the condition number is less than machine
76 *> precision, INFO = N+1 is returned as a warning, but the routine
77 *> still goes on to solve for X and compute error bounds as
80 *> 4. The system of equations is solved for X using the factored form
83 *> 5. Iterative refinement is applied to improve the computed solution
84 *> matrix and calculate error bounds and backward error estimates
87 *> 6. If equilibration was used, the matrix X is premultiplied by
88 *> diag(S) so that it solves the original system before
97 *> FACT is CHARACTER*1
98 *> Specifies whether or not the factored form of the matrix A is
99 *> supplied on entry, and if not, whether the matrix A should be
100 *> equilibrated before it is factored.
101 *> = 'F': On entry, AFP contains the factored form of A.
102 *> If EQUED = 'Y', the matrix A has been equilibrated
103 *> with scaling factors given by S. AP and AFP will not
105 *> = 'N': The matrix A will be copied to AFP and factored.
106 *> = 'E': The matrix A will be equilibrated if necessary, then
107 *> copied to AFP and factored.
112 *> UPLO is CHARACTER*1
113 *> = 'U': Upper triangle of A is stored;
114 *> = 'L': Lower triangle of A is stored.
120 *> The number of linear equations, i.e., the order of the
127 *> The number of right hand sides, i.e., the number of columns
128 *> of the matrices B and X. NRHS >= 0.
133 *> AP is DOUBLE PRECISION array, dimension (N*(N+1)/2)
134 *> On entry, the upper or lower triangle of the symmetric matrix
135 *> A, packed columnwise in a linear array, except if FACT = 'F'
136 *> and EQUED = 'Y', then A must contain the equilibrated matrix
137 *> diag(S)*A*diag(S). The j-th column of A is stored in the
138 *> array AP as follows:
139 *> if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
140 *> if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.
141 *> See below for further details. A is not modified if
142 *> FACT = 'F' or 'N', or if FACT = 'E' and EQUED = 'N' on exit.
144 *> On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
145 *> diag(S)*A*diag(S).
148 *> \param[in,out] AFP
150 *> AFP is DOUBLE PRECISION array, dimension
152 *> If FACT = 'F', then AFP is an input argument and on entry
153 *> contains the triangular factor U or L from the Cholesky
154 *> factorization A = U**T*U or A = L*L**T, in the same storage
155 *> format as A. If EQUED .ne. 'N', then AFP is the factored
156 *> form of the equilibrated matrix A.
158 *> If FACT = 'N', then AFP is an output argument and on exit
159 *> returns the triangular factor U or L from the Cholesky
160 *> factorization A = U**T * U or A = L * L**T of the original
163 *> If FACT = 'E', then AFP is an output argument and on exit
164 *> returns the triangular factor U or L from the Cholesky
165 *> factorization A = U**T * U or A = L * L**T of the equilibrated
166 *> matrix A (see the description of AP for the form of the
167 *> equilibrated matrix).
170 *> \param[in,out] EQUED
172 *> EQUED is CHARACTER*1
173 *> Specifies the form of equilibration that was done.
174 *> = 'N': No equilibration (always true if FACT = 'N').
175 *> = 'Y': Equilibration was done, i.e., A has been replaced by
176 *> diag(S) * A * diag(S).
177 *> EQUED is an input argument if FACT = 'F'; otherwise, it is an
183 *> S is DOUBLE PRECISION array, dimension (N)
184 *> The scale factors for A; not accessed if EQUED = 'N'. S is
185 *> an input argument if FACT = 'F'; otherwise, S is an output
186 *> argument. If FACT = 'F' and EQUED = 'Y', each element of S
192 *> B is DOUBLE PRECISION array, dimension (LDB,NRHS)
193 *> On entry, the N-by-NRHS right hand side matrix B.
194 *> On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y',
195 *> B is overwritten by diag(S) * B.
201 *> The leading dimension of the array B. LDB >= max(1,N).
206 *> X is DOUBLE PRECISION array, dimension (LDX,NRHS)
207 *> If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X to
208 *> the original system of equations. Note that if EQUED = 'Y',
209 *> A and B are modified on exit, and the solution to the
210 *> equilibrated system is inv(diag(S))*X.
216 *> The leading dimension of the array X. LDX >= max(1,N).
221 *> RCOND is DOUBLE PRECISION
222 *> The estimate of the reciprocal condition number of the matrix
223 *> A after equilibration (if done). If RCOND is less than the
224 *> machine precision (in particular, if RCOND = 0), the matrix
225 *> is singular to working precision. This condition is
226 *> indicated by a return code of INFO > 0.
231 *> FERR is DOUBLE PRECISION array, dimension (NRHS)
232 *> The estimated forward error bound for each solution vector
233 *> X(j) (the j-th column of the solution matrix X).
234 *> If XTRUE is the true solution corresponding to X(j), FERR(j)
235 *> is an estimated upper bound for the magnitude of the largest
236 *> element in (X(j) - XTRUE) divided by the magnitude of the
237 *> largest element in X(j). The estimate is as reliable as
238 *> the estimate for RCOND, and is almost always a slight
239 *> overestimate of the true error.
244 *> BERR is DOUBLE PRECISION array, dimension (NRHS)
245 *> The componentwise relative backward error of each solution
246 *> vector X(j) (i.e., the smallest relative change in
247 *> any element of A or B that makes X(j) an exact solution).
252 *> WORK is DOUBLE PRECISION array, dimension (3*N)
257 *> IWORK is INTEGER array, dimension (N)
263 *> = 0: successful exit
264 *> < 0: if INFO = -i, the i-th argument had an illegal value
265 *> > 0: if INFO = i, and i is
266 *> <= N: the leading minor of order i of A is
267 *> not positive definite, so the factorization
268 *> could not be completed, and the solution has not
269 *> been computed. RCOND = 0 is returned.
270 *> = N+1: U is nonsingular, but RCOND is less than machine
271 *> precision, meaning that the matrix is singular
272 *> to working precision. Nevertheless, the
273 *> solution and error bounds are computed because
274 *> there are a number of situations where the
275 *> computed solution can be more accurate than the
276 *> value of RCOND would suggest.
282 *> \author Univ. of Tennessee
283 *> \author Univ. of California Berkeley
284 *> \author Univ. of Colorado Denver
289 *> \ingroup doubleOTHERsolve
291 *> \par Further Details:
292 * =====================
296 *> The packed storage scheme is illustrated by the following example
297 *> when N = 4, UPLO = 'U':
299 *> Two-dimensional storage of the symmetric matrix A:
303 *> a33 a34 (aij = conjg(aji))
306 *> Packed storage of the upper triangle of A:
308 *> AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]
311 * =====================================================================
312 SUBROUTINE DPPSVX( FACT, UPLO, N, NRHS, AP, AFP, EQUED, S, B, LDB,
313 $ X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO )
315 * -- LAPACK driver routine (version 3.4.1) --
316 * -- LAPACK is a software package provided by Univ. of Tennessee, --
317 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
320 * .. Scalar Arguments ..
321 CHARACTER EQUED, FACT, UPLO
322 INTEGER INFO, LDB, LDX, N, NRHS
323 DOUBLE PRECISION RCOND
325 * .. Array Arguments ..
327 DOUBLE PRECISION AFP( * ), AP( * ), B( LDB, * ), BERR( * ),
328 $ FERR( * ), S( * ), WORK( * ), X( LDX, * )
331 * =====================================================================
334 DOUBLE PRECISION ZERO, ONE
335 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
337 * .. Local Scalars ..
338 LOGICAL EQUIL, NOFACT, RCEQU
340 DOUBLE PRECISION AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM
342 * .. External Functions ..
344 DOUBLE PRECISION DLAMCH, DLANSP
345 EXTERNAL LSAME, DLAMCH, DLANSP
347 * .. External Subroutines ..
348 EXTERNAL DCOPY, DLACPY, DLAQSP, DPPCON, DPPEQU, DPPRFS,
349 $ DPPTRF, DPPTRS, XERBLA
351 * .. Intrinsic Functions ..
354 * .. Executable Statements ..
357 NOFACT = LSAME( FACT, 'N' )
358 EQUIL = LSAME( FACT, 'E' )
359 IF( NOFACT .OR. EQUIL ) THEN
363 RCEQU = LSAME( EQUED, 'Y' )
364 SMLNUM = DLAMCH( 'Safe minimum' )
365 BIGNUM = ONE / SMLNUM
368 * Test the input parameters.
370 IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.LSAME( FACT, 'F' ) )
373 ELSE IF( .NOT.LSAME( UPLO, 'U' ) .AND. .NOT.LSAME( UPLO, 'L' ) )
376 ELSE IF( N.LT.0 ) THEN
378 ELSE IF( NRHS.LT.0 ) THEN
380 ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
381 $ ( RCEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
388 SMIN = MIN( SMIN, S( J ) )
389 SMAX = MAX( SMAX, S( J ) )
391 IF( SMIN.LE.ZERO ) THEN
393 ELSE IF( N.GT.0 ) THEN
394 SCOND = MAX( SMIN, SMLNUM ) / MIN( SMAX, BIGNUM )
400 IF( LDB.LT.MAX( 1, N ) ) THEN
402 ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
409 CALL XERBLA( 'DPPSVX', -INFO )
415 * Compute row and column scalings to equilibrate the matrix A.
417 CALL DPPEQU( UPLO, N, AP, S, SCOND, AMAX, INFEQU )
418 IF( INFEQU.EQ.0 ) THEN
420 * Equilibrate the matrix.
422 CALL DLAQSP( UPLO, N, AP, S, SCOND, AMAX, EQUED )
423 RCEQU = LSAME( EQUED, 'Y' )
427 * Scale the right-hand side.
432 B( I, J ) = S( I )*B( I, J )
437 IF( NOFACT .OR. EQUIL ) THEN
439 * Compute the Cholesky factorization A = U**T * U or A = L * L**T.
441 CALL DCOPY( N*( N+1 ) / 2, AP, 1, AFP, 1 )
442 CALL DPPTRF( UPLO, N, AFP, INFO )
444 * Return if INFO is non-zero.
452 * Compute the norm of the matrix A.
454 ANORM = DLANSP( 'I', UPLO, N, AP, WORK )
456 * Compute the reciprocal of the condition number of A.
458 CALL DPPCON( UPLO, N, AFP, ANORM, RCOND, WORK, IWORK, INFO )
460 * Compute the solution matrix X.
462 CALL DLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
463 CALL DPPTRS( UPLO, N, NRHS, AFP, X, LDX, INFO )
465 * Use iterative refinement to improve the computed solution and
466 * compute error bounds and backward error estimates for it.
468 CALL DPPRFS( UPLO, N, NRHS, AP, AFP, B, LDB, X, LDX, FERR, BERR,
469 $ WORK, IWORK, INFO )
471 * Transform the solution matrix X to a solution of the original
477 X( I, J ) = S( I )*X( I, J )
481 FERR( J ) = FERR( J ) / SCOND
485 * Set INFO = N+1 if the matrix is singular to working precision.
487 IF( RCOND.LT.DLAMCH( 'Epsilon' ) )