1 *> \brief <b> SPOSVXX computes the solution to system of linear equations A * X = B for PO matrices</b>
3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download SPOSVXX + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sposvxx.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sposvxx.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sposvxx.f">
21 * SUBROUTINE SPOSVXX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED,
22 * S, B, LDB, X, LDX, RCOND, RPVGRW, BERR,
23 * N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP,
24 * NPARAMS, PARAMS, WORK, IWORK, INFO )
26 * .. Scalar Arguments ..
27 * CHARACTER EQUED, FACT, UPLO
28 * INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
32 * .. Array Arguments ..
34 * REAL A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
35 * $ X( LDX, * ), WORK( * )
36 * REAL S( * ), PARAMS( * ), BERR( * ),
37 * $ ERR_BNDS_NORM( NRHS, * ),
38 * $ ERR_BNDS_COMP( NRHS, * )
47 *> SPOSVXX uses the Cholesky factorization A = U**T*U or A = L*L**T
48 *> to compute the solution to a real system of linear equations
49 *> A * X = B, where A is an N-by-N symmetric positive definite matrix
50 *> and X and B are N-by-NRHS matrices.
52 *> If requested, both normwise and maximum componentwise error bounds
53 *> are returned. SPOSVXX will return a solution with a tiny
54 *> guaranteed error (O(eps) where eps is the working machine
55 *> precision) unless the matrix is very ill-conditioned, in which
56 *> case a warning is returned. Relevant condition numbers also are
57 *> calculated and returned.
59 *> SPOSVXX accepts user-provided factorizations and equilibration
60 *> factors; see the definitions of the FACT and EQUED options.
61 *> Solving with refinement and using a factorization from a previous
62 *> SPOSVXX call will also produce a solution with either O(eps)
63 *> errors or warnings, but we cannot make that claim for general
64 *> user-provided factorizations and equilibration factors if they
65 *> differ from what SPOSVXX would itself produce.
73 *> The following steps are performed:
75 *> 1. If FACT = 'E', real scaling factors are computed to equilibrate
78 *> diag(S)*A*diag(S) *inv(diag(S))*X = diag(S)*B
80 *> Whether or not the system will be equilibrated depends on the
81 *> scaling of the matrix A, but if equilibration is used, A is
82 *> overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
84 *> 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
85 *> factor the matrix A (after equilibration if FACT = 'E') as
86 *> A = U**T* U, if UPLO = 'U', or
87 *> A = L * L**T, if UPLO = 'L',
88 *> where U is an upper triangular matrix and L is a lower triangular
91 *> 3. If the leading i-by-i principal minor is not positive definite,
92 *> then the routine returns with INFO = i. Otherwise, the factored
93 *> form of A is used to estimate the condition number of the matrix
94 *> A (see argument RCOND). If the reciprocal of the condition number
95 *> is less than machine precision, the routine still goes on to solve
96 *> for X and compute error bounds as described below.
98 *> 4. The system of equations is solved for X using the factored form
101 *> 5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
102 *> the routine will use iterative refinement to try to get a small
103 *> error and error bounds. Refinement calculates the residual to at
104 *> least twice the working precision.
106 *> 6. If equilibration was used, the matrix X is premultiplied by
107 *> diag(S) so that it solves the original system before
115 *> Some optional parameters are bundled in the PARAMS array. These
116 *> settings determine how refinement is performed, but often the
117 *> defaults are acceptable. If the defaults are acceptable, users
118 *> can pass NPARAMS = 0 which prevents the source code from accessing
119 *> the PARAMS argument.
124 *> FACT is CHARACTER*1
125 *> Specifies whether or not the factored form of the matrix A is
126 *> supplied on entry, and if not, whether the matrix A should be
127 *> equilibrated before it is factored.
128 *> = 'F': On entry, AF contains the factored form of A.
129 *> If EQUED is not 'N', the matrix A has been
130 *> equilibrated with scaling factors given by S.
131 *> A and AF are not modified.
132 *> = 'N': The matrix A will be copied to AF and factored.
133 *> = 'E': The matrix A will be equilibrated if necessary, then
134 *> copied to AF and factored.
139 *> UPLO is CHARACTER*1
140 *> = 'U': Upper triangle of A is stored;
141 *> = 'L': Lower triangle of A is stored.
147 *> The number of linear equations, i.e., the order of the
154 *> The number of right hand sides, i.e., the number of columns
155 *> of the matrices B and X. NRHS >= 0.
160 *> A is REAL array, dimension (LDA,N)
161 *> On entry, the symmetric matrix A, except if FACT = 'F' and EQUED =
162 *> 'Y', then A must contain the equilibrated matrix
163 *> diag(S)*A*diag(S). If UPLO = 'U', the leading N-by-N upper
164 *> triangular part of A contains the upper triangular part of the
165 *> matrix A, and the strictly lower triangular part of A is not
166 *> referenced. If UPLO = 'L', the leading N-by-N lower triangular
167 *> part of A contains the lower triangular part of the matrix A, and
168 *> the strictly upper triangular part of A is not referenced. A is
169 *> not modified if FACT = 'F' or 'N', or if FACT = 'E' and EQUED =
172 *> On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
173 *> diag(S)*A*diag(S).
179 *> The leading dimension of the array A. LDA >= max(1,N).
184 *> AF is REAL array, dimension (LDAF,N)
185 *> If FACT = 'F', then AF is an input argument and on entry
186 *> contains the triangular factor U or L from the Cholesky
187 *> factorization A = U**T*U or A = L*L**T, in the same storage
188 *> format as A. If EQUED .ne. 'N', then AF is the factored
189 *> form of the equilibrated matrix diag(S)*A*diag(S).
191 *> If FACT = 'N', then AF is an output argument and on exit
192 *> returns the triangular factor U or L from the Cholesky
193 *> factorization A = U**T*U or A = L*L**T of the original
196 *> If FACT = 'E', then AF is an output argument and on exit
197 *> returns the triangular factor U or L from the Cholesky
198 *> factorization A = U**T*U or A = L*L**T of the equilibrated
199 *> matrix A (see the description of A for the form of the
200 *> equilibrated matrix).
206 *> The leading dimension of the array AF. LDAF >= max(1,N).
209 *> \param[in,out] EQUED
211 *> EQUED is CHARACTER*1
212 *> Specifies the form of equilibration that was done.
213 *> = 'N': No equilibration (always true if FACT = 'N').
214 *> = 'Y': Both row and column equilibration, i.e., A has been
215 *> replaced by diag(S) * A * diag(S).
216 *> EQUED is an input argument if FACT = 'F'; otherwise, it is an
222 *> S is REAL array, dimension (N)
223 *> The row scale factors for A. If EQUED = 'Y', A is multiplied on
224 *> the left and right by diag(S). S is an input argument if FACT =
225 *> 'F'; otherwise, S is an output argument. If FACT = 'F' and EQUED
226 *> = 'Y', each element of S must be positive. If S is output, each
227 *> element of S is a power of the radix. If S is input, each element
228 *> of S should be a power of the radix to ensure a reliable solution
229 *> and error estimates. Scaling by powers of the radix does not cause
230 *> rounding errors unless the result underflows or overflows.
231 *> Rounding errors during scaling lead to refining with a matrix that
232 *> is not equivalent to the input matrix, producing error estimates
233 *> that may not be reliable.
238 *> B is REAL array, dimension (LDB,NRHS)
239 *> On entry, the N-by-NRHS right hand side matrix B.
241 *> if EQUED = 'N', B is not modified;
242 *> if EQUED = 'Y', B is overwritten by diag(S)*B;
248 *> The leading dimension of the array B. LDB >= max(1,N).
253 *> X is REAL array, dimension (LDX,NRHS)
254 *> If INFO = 0, the N-by-NRHS solution matrix X to the original
255 *> system of equations. Note that A and B are modified on exit if
256 *> EQUED .ne. 'N', and the solution to the equilibrated system is
263 *> The leading dimension of the array X. LDX >= max(1,N).
269 *> Reciprocal scaled condition number. This is an estimate of the
270 *> reciprocal Skeel condition number of the matrix A after
271 *> equilibration (if done). If this is less than the machine
272 *> precision (in particular, if it is zero), the matrix is singular
273 *> to working precision. Note that the error may still be small even
274 *> if this number is very small and the matrix appears ill-
278 *> \param[out] RPVGRW
281 *> Reciprocal pivot growth. On exit, this contains the reciprocal
282 *> pivot growth factor norm(A)/norm(U). The "max absolute element"
283 *> norm is used. If this is much less than 1, then the stability of
284 *> the LU factorization of the (equilibrated) matrix A could be poor.
285 *> This also means that the solution X, estimated condition numbers,
286 *> and error bounds could be unreliable. If factorization fails with
287 *> 0<INFO<=N, then this contains the reciprocal pivot growth factor
288 *> for the leading INFO columns of A.
293 *> BERR is REAL array, dimension (NRHS)
294 *> Componentwise relative backward error. This is the
295 *> componentwise relative backward error of each solution vector X(j)
296 *> (i.e., the smallest relative change in any element of A or B that
297 *> makes X(j) an exact solution).
300 *> \param[in] N_ERR_BNDS
302 *> N_ERR_BNDS is INTEGER
303 *> Number of error bounds to return for each right hand side
304 *> and each type (normwise or componentwise). See ERR_BNDS_NORM and
305 *> ERR_BNDS_COMP below.
308 *> \param[out] ERR_BNDS_NORM
310 *> ERR_BNDS_NORM is REAL array, dimension (NRHS, N_ERR_BNDS)
311 *> For each right-hand side, this array contains information about
312 *> various error bounds and condition numbers corresponding to the
313 *> normwise relative error, which is defined as follows:
315 *> Normwise relative error in the ith solution vector:
316 *> max_j (abs(XTRUE(j,i) - X(j,i)))
317 *> ------------------------------
320 *> The array is indexed by the type of error information as described
321 *> below. There currently are up to three pieces of information
324 *> The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
327 *> The second index in ERR_BNDS_NORM(:,err) contains the following
329 *> err = 1 "Trust/don't trust" boolean. Trust the answer if the
330 *> reciprocal condition number is less than the threshold
331 *> sqrt(n) * slamch('Epsilon').
333 *> err = 2 "Guaranteed" error bound: The estimated forward error,
334 *> almost certainly within a factor of 10 of the true error
335 *> so long as the next entry is greater than the threshold
336 *> sqrt(n) * slamch('Epsilon'). This error bound should only
337 *> be trusted if the previous boolean is true.
339 *> err = 3 Reciprocal condition number: Estimated normwise
340 *> reciprocal condition number. Compared with the threshold
341 *> sqrt(n) * slamch('Epsilon') to determine if the error
342 *> estimate is "guaranteed". These reciprocal condition
343 *> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
344 *> appropriately scaled matrix Z.
345 *> Let Z = S*A, where S scales each row by a power of the
346 *> radix so all absolute row sums of Z are approximately 1.
348 *> See Lapack Working Note 165 for further details and extra
352 *> \param[out] ERR_BNDS_COMP
354 *> ERR_BNDS_COMP is REAL array, dimension (NRHS, N_ERR_BNDS)
355 *> For each right-hand side, this array contains information about
356 *> various error bounds and condition numbers corresponding to the
357 *> componentwise relative error, which is defined as follows:
359 *> Componentwise relative error in the ith solution vector:
360 *> abs(XTRUE(j,i) - X(j,i))
361 *> max_j ----------------------
364 *> The array is indexed by the right-hand side i (on which the
365 *> componentwise relative error depends), and the type of error
366 *> information as described below. There currently are up to three
367 *> pieces of information returned for each right-hand side. If
368 *> componentwise accuracy is not requested (PARAMS(3) = 0.0), then
369 *> ERR_BNDS_COMP is not accessed. If N_ERR_BNDS .LT. 3, then at most
370 *> the first (:,N_ERR_BNDS) entries are returned.
372 *> The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
375 *> The second index in ERR_BNDS_COMP(:,err) contains the following
377 *> err = 1 "Trust/don't trust" boolean. Trust the answer if the
378 *> reciprocal condition number is less than the threshold
379 *> sqrt(n) * slamch('Epsilon').
381 *> err = 2 "Guaranteed" error bound: The estimated forward error,
382 *> almost certainly within a factor of 10 of the true error
383 *> so long as the next entry is greater than the threshold
384 *> sqrt(n) * slamch('Epsilon'). This error bound should only
385 *> be trusted if the previous boolean is true.
387 *> err = 3 Reciprocal condition number: Estimated componentwise
388 *> reciprocal condition number. Compared with the threshold
389 *> sqrt(n) * slamch('Epsilon') to determine if the error
390 *> estimate is "guaranteed". These reciprocal condition
391 *> numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
392 *> appropriately scaled matrix Z.
393 *> Let Z = S*(A*diag(x)), where x is the solution for the
394 *> current right-hand side and S scales each row of
395 *> A*diag(x) by a power of the radix so all absolute row
396 *> sums of Z are approximately 1.
398 *> See Lapack Working Note 165 for further details and extra
402 *> \param[in] NPARAMS
404 *> NPARAMS is INTEGER
405 *> Specifies the number of parameters set in PARAMS. If .LE. 0, the
406 *> PARAMS array is never referenced and default values are used.
409 *> \param[in,out] PARAMS
411 *> PARAMS is REAL array, dimension NPARAMS
412 *> Specifies algorithm parameters. If an entry is .LT. 0.0, then
413 *> that entry will be filled with default value used for that
414 *> parameter. Only positions up to NPARAMS are accessed; defaults
415 *> are used for higher-numbered parameters.
417 *> PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
418 *> refinement or not.
420 *> = 0.0 : No refinement is performed, and no error bounds are
422 *> = 1.0 : Use the double-precision refinement algorithm,
423 *> possibly with doubled-single computations if the
424 *> compilation environment does not support DOUBLE
426 *> (other values are reserved for future use)
428 *> PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
429 *> computations allowed for refinement.
431 *> Aggressive: Set to 100 to permit convergence using approximate
432 *> factorizations or factorizations other than LU. If
433 *> the factorization uses a technique other than
434 *> Gaussian elimination, the guarantees in
435 *> err_bnds_norm and err_bnds_comp may no longer be
438 *> PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
439 *> will attempt to find a solution with small componentwise
440 *> relative error in the double-precision algorithm. Positive
441 *> is true, 0.0 is false.
442 *> Default: 1.0 (attempt componentwise convergence)
447 *> WORK is REAL array, dimension (4*N)
452 *> IWORK is INTEGER array, dimension (N)
458 *> = 0: Successful exit. The solution to every right-hand side is
460 *> < 0: If INFO = -i, the i-th argument had an illegal value
461 *> > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization
462 *> has been completed, but the factor U is exactly singular, so
463 *> the solution and error bounds could not be computed. RCOND = 0
465 *> = N+J: The solution corresponding to the Jth right-hand side is
466 *> not guaranteed. The solutions corresponding to other right-
467 *> hand sides K with K > J may not be guaranteed as well, but
468 *> only the first such right-hand side is reported. If a small
469 *> componentwise error is not requested (PARAMS(3) = 0.0) then
470 *> the Jth right-hand side is the first with a normwise error
471 *> bound that is not guaranteed (the smallest J such
472 *> that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
473 *> the Jth right-hand side is the first with either a normwise or
474 *> componentwise error bound that is not guaranteed (the smallest
475 *> J such that either ERR_BNDS_NORM(J,1) = 0.0 or
476 *> ERR_BNDS_COMP(J,1) = 0.0). See the definition of
477 *> ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
478 *> about all of the right-hand sides check ERR_BNDS_NORM or
485 *> \author Univ. of Tennessee
486 *> \author Univ. of California Berkeley
487 *> \author Univ. of Colorado Denver
492 *> \ingroup realPOsolve
494 * =====================================================================
495 SUBROUTINE SPOSVXX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED,
496 $ S, B, LDB, X, LDX, RCOND, RPVGRW, BERR,
497 $ N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP,
498 $ NPARAMS, PARAMS, WORK, IWORK, INFO )
500 * -- LAPACK driver routine (version 3.4.1) --
501 * -- LAPACK is a software package provided by Univ. of Tennessee, --
502 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
505 * .. Scalar Arguments ..
506 CHARACTER EQUED, FACT, UPLO
507 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
511 * .. Array Arguments ..
513 REAL A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
514 $ X( LDX, * ), WORK( * )
515 REAL S( * ), PARAMS( * ), BERR( * ),
516 $ ERR_BNDS_NORM( NRHS, * ),
517 $ ERR_BNDS_COMP( NRHS, * )
520 * ==================================================================
524 PARAMETER ( ZERO = 0.0E+0, ONE = 1.0E+0 )
525 INTEGER FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
526 INTEGER RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
527 INTEGER CMP_ERR_I, PIV_GROWTH_I
528 PARAMETER ( FINAL_NRM_ERR_I = 1, FINAL_CMP_ERR_I = 2,
530 PARAMETER ( RCOND_I = 4, NRM_RCOND_I = 5, NRM_ERR_I = 6 )
531 PARAMETER ( CMP_RCOND_I = 7, CMP_ERR_I = 8,
534 * .. Local Scalars ..
535 LOGICAL EQUIL, NOFACT, RCEQU
537 REAL AMAX, BIGNUM, SMIN, SMAX,
540 * .. External Functions ..
541 EXTERNAL LSAME, SLAMCH, SLA_PORPVGRW
543 REAL SLAMCH, SLA_PORPVGRW
545 * .. External Subroutines ..
546 EXTERNAL SPOEQUB, SPOTRF, SPOTRS, SLACPY, SLAQSY,
547 $ XERBLA, SLASCL2, SPORFSX
549 * .. Intrinsic Functions ..
552 * .. Executable Statements ..
555 NOFACT = LSAME( FACT, 'N' )
556 EQUIL = LSAME( FACT, 'E' )
557 SMLNUM = SLAMCH( 'Safe minimum' )
558 BIGNUM = ONE / SMLNUM
559 IF( NOFACT .OR. EQUIL ) THEN
563 RCEQU = LSAME( EQUED, 'Y' )
566 * Default is failure. If an input parameter is wrong or
567 * factorization fails, make everything look horrible. Only the
568 * pivot growth is set here, the rest is initialized in SPORFSX.
572 * Test the input parameters. PARAMS is not tested until SPORFSX.
574 IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.
575 $ LSAME( FACT, 'F' ) ) THEN
577 ELSE IF( .NOT.LSAME( UPLO, 'U' ) .AND.
578 $ .NOT.LSAME( UPLO, 'L' ) ) THEN
580 ELSE IF( N.LT.0 ) THEN
582 ELSE IF( NRHS.LT.0 ) THEN
584 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
586 ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN
588 ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
589 $ ( RCEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
596 SMIN = MIN( SMIN, S( J ) )
597 SMAX = MAX( SMAX, S( J ) )
599 IF( SMIN.LE.ZERO ) THEN
601 ELSE IF( N.GT.0 ) THEN
602 SCOND = MAX( SMIN, SMLNUM ) / MIN( SMAX, BIGNUM )
608 IF( LDB.LT.MAX( 1, N ) ) THEN
610 ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
617 CALL XERBLA( 'SPOSVXX', -INFO )
623 * Compute row and column scalings to equilibrate the matrix A.
625 CALL SPOEQUB( N, A, LDA, S, SCOND, AMAX, INFEQU )
626 IF( INFEQU.EQ.0 ) THEN
628 * Equilibrate the matrix.
630 CALL SLAQSY( UPLO, N, A, LDA, S, SCOND, AMAX, EQUED )
631 RCEQU = LSAME( EQUED, 'Y' )
635 * Scale the right-hand side.
637 IF( RCEQU ) CALL SLASCL2( N, NRHS, S, B, LDB )
639 IF( NOFACT .OR. EQUIL ) THEN
641 * Compute the Cholesky factorization of A.
643 CALL SLACPY( UPLO, N, N, A, LDA, AF, LDAF )
644 CALL SPOTRF( UPLO, N, AF, LDAF, INFO )
646 * Return if INFO is non-zero.
650 * Pivot in column INFO is exactly 0
651 * Compute the reciprocal pivot growth factor of the
652 * leading rank-deficient INFO columns of A.
654 RPVGRW = SLA_PORPVGRW( UPLO, INFO, A, LDA, AF, LDAF, WORK )
659 * Compute the reciprocal growth factor RPVGRW.
661 RPVGRW = SLA_PORPVGRW( UPLO, N, A, LDA, AF, LDAF, WORK )
663 * Compute the solution matrix X.
665 CALL SLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
666 CALL SPOTRS( UPLO, N, NRHS, AF, LDAF, X, LDX, INFO )
668 * Use iterative refinement to improve the computed solution and
669 * compute error bounds and backward error estimates for it.
671 CALL SPORFSX( UPLO, EQUED, N, NRHS, A, LDA, AF, LDAF,
672 $ S, B, LDB, X, LDX, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_NORM,
673 $ ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, IWORK, INFO )
679 CALL SLASCL2 ( N, NRHS, S, X, LDX )