LAPACK 3.3.1
Linear Algebra PACKage
|
00001 SUBROUTINE SPOSVX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED, 00002 $ S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, 00003 $ IWORK, INFO ) 00004 * 00005 * -- LAPACK driver routine (version 3.3.1) -- 00006 * -- LAPACK is a software package provided by Univ. of Tennessee, -- 00007 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 00008 * -- April 2011 -- 00009 * 00010 * .. Scalar Arguments .. 00011 CHARACTER EQUED, FACT, UPLO 00012 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS 00013 REAL RCOND 00014 * .. 00015 * .. Array Arguments .. 00016 INTEGER IWORK( * ) 00017 REAL A( LDA, * ), AF( LDAF, * ), B( LDB, * ), 00018 $ BERR( * ), FERR( * ), S( * ), WORK( * ), 00019 $ X( LDX, * ) 00020 * .. 00021 * 00022 * Purpose 00023 * ======= 00024 * 00025 * SPOSVX uses the Cholesky factorization A = U**T*U or A = L*L**T to 00026 * compute the solution to a real system of linear equations 00027 * A * X = B, 00028 * where A is an N-by-N symmetric positive definite matrix and X and B 00029 * are N-by-NRHS matrices. 00030 * 00031 * Error bounds on the solution and a condition estimate are also 00032 * provided. 00033 * 00034 * Description 00035 * =========== 00036 * 00037 * The following steps are performed: 00038 * 00039 * 1. If FACT = 'E', real scaling factors are computed to equilibrate 00040 * the system: 00041 * diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * B 00042 * Whether or not the system will be equilibrated depends on the 00043 * scaling of the matrix A, but if equilibration is used, A is 00044 * overwritten by diag(S)*A*diag(S) and B by diag(S)*B. 00045 * 00046 * 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to 00047 * factor the matrix A (after equilibration if FACT = 'E') as 00048 * A = U**T* U, if UPLO = 'U', or 00049 * A = L * L**T, if UPLO = 'L', 00050 * where U is an upper triangular matrix and L is a lower triangular 00051 * matrix. 00052 * 00053 * 3. If the leading i-by-i principal minor is not positive definite, 00054 * then the routine returns with INFO = i. Otherwise, the factored 00055 * form of A is used to estimate the condition number of the matrix 00056 * A. If the reciprocal of the condition number is less than machine 00057 * precision, INFO = N+1 is returned as a warning, but the routine 00058 * still goes on to solve for X and compute error bounds as 00059 * described below. 00060 * 00061 * 4. The system of equations is solved for X using the factored form 00062 * of A. 00063 * 00064 * 5. Iterative refinement is applied to improve the computed solution 00065 * matrix and calculate error bounds and backward error estimates 00066 * for it. 00067 * 00068 * 6. If equilibration was used, the matrix X is premultiplied by 00069 * diag(S) so that it solves the original system before 00070 * equilibration. 00071 * 00072 * Arguments 00073 * ========= 00074 * 00075 * FACT (input) CHARACTER*1 00076 * Specifies whether or not the factored form of the matrix A is 00077 * supplied on entry, and if not, whether the matrix A should be 00078 * equilibrated before it is factored. 00079 * = 'F': On entry, AF contains the factored form of A. 00080 * If EQUED = 'Y', the matrix A has been equilibrated 00081 * with scaling factors given by S. A and AF will not 00082 * be modified. 00083 * = 'N': The matrix A will be copied to AF and factored. 00084 * = 'E': The matrix A will be equilibrated if necessary, then 00085 * copied to AF and factored. 00086 * 00087 * UPLO (input) CHARACTER*1 00088 * = 'U': Upper triangle of A is stored; 00089 * = 'L': Lower triangle of A is stored. 00090 * 00091 * N (input) INTEGER 00092 * The number of linear equations, i.e., the order of the 00093 * matrix A. N >= 0. 00094 * 00095 * NRHS (input) INTEGER 00096 * The number of right hand sides, i.e., the number of columns 00097 * of the matrices B and X. NRHS >= 0. 00098 * 00099 * A (input/output) REAL array, dimension (LDA,N) 00100 * On entry, the symmetric matrix A, except if FACT = 'F' and 00101 * EQUED = 'Y', then A must contain the equilibrated matrix 00102 * diag(S)*A*diag(S). If UPLO = 'U', the leading 00103 * N-by-N upper triangular part of A contains the upper 00104 * triangular part of the matrix A, and the strictly lower 00105 * triangular part of A is not referenced. If UPLO = 'L', the 00106 * leading N-by-N lower triangular part of A contains the lower 00107 * triangular part of the matrix A, and the strictly upper 00108 * triangular part of A is not referenced. A is not modified if 00109 * FACT = 'F' or 'N', or if FACT = 'E' and EQUED = 'N' on exit. 00110 * 00111 * On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by 00112 * diag(S)*A*diag(S). 00113 * 00114 * LDA (input) INTEGER 00115 * The leading dimension of the array A. LDA >= max(1,N). 00116 * 00117 * AF (input or output) REAL array, dimension (LDAF,N) 00118 * If FACT = 'F', then AF is an input argument and on entry 00119 * contains the triangular factor U or L from the Cholesky 00120 * factorization A = U**T*U or A = L*L**T, in the same storage 00121 * format as A. If EQUED .ne. 'N', then AF is the factored form 00122 * of the equilibrated matrix diag(S)*A*diag(S). 00123 * 00124 * If FACT = 'N', then AF is an output argument and on exit 00125 * returns the triangular factor U or L from the Cholesky 00126 * factorization A = U**T*U or A = L*L**T of the original 00127 * matrix A. 00128 * 00129 * If FACT = 'E', then AF is an output argument and on exit 00130 * returns the triangular factor U or L from the Cholesky 00131 * factorization A = U**T*U or A = L*L**T of the equilibrated 00132 * matrix A (see the description of A for the form of the 00133 * equilibrated matrix). 00134 * 00135 * LDAF (input) INTEGER 00136 * The leading dimension of the array AF. LDAF >= max(1,N). 00137 * 00138 * EQUED (input or output) CHARACTER*1 00139 * Specifies the form of equilibration that was done. 00140 * = 'N': No equilibration (always true if FACT = 'N'). 00141 * = 'Y': Equilibration was done, i.e., A has been replaced by 00142 * diag(S) * A * diag(S). 00143 * EQUED is an input argument if FACT = 'F'; otherwise, it is an 00144 * output argument. 00145 * 00146 * S (input or output) REAL array, dimension (N) 00147 * The scale factors for A; not accessed if EQUED = 'N'. S is 00148 * an input argument if FACT = 'F'; otherwise, S is an output 00149 * argument. If FACT = 'F' and EQUED = 'Y', each element of S 00150 * must be positive. 00151 * 00152 * B (input/output) REAL array, dimension (LDB,NRHS) 00153 * On entry, the N-by-NRHS right hand side matrix B. 00154 * On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y', 00155 * B is overwritten by diag(S) * B. 00156 * 00157 * LDB (input) INTEGER 00158 * The leading dimension of the array B. LDB >= max(1,N). 00159 * 00160 * X (output) REAL array, dimension (LDX,NRHS) 00161 * If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X to 00162 * the original system of equations. Note that if EQUED = 'Y', 00163 * A and B are modified on exit, and the solution to the 00164 * equilibrated system is inv(diag(S))*X. 00165 * 00166 * LDX (input) INTEGER 00167 * The leading dimension of the array X. LDX >= max(1,N). 00168 * 00169 * RCOND (output) REAL 00170 * The estimate of the reciprocal condition number of the matrix 00171 * A after equilibration (if done). If RCOND is less than the 00172 * machine precision (in particular, if RCOND = 0), the matrix 00173 * is singular to working precision. This condition is 00174 * indicated by a return code of INFO > 0. 00175 * 00176 * FERR (output) REAL array, dimension (NRHS) 00177 * The estimated forward error bound for each solution vector 00178 * X(j) (the j-th column of the solution matrix X). 00179 * If XTRUE is the true solution corresponding to X(j), FERR(j) 00180 * is an estimated upper bound for the magnitude of the largest 00181 * element in (X(j) - XTRUE) divided by the magnitude of the 00182 * largest element in X(j). The estimate is as reliable as 00183 * the estimate for RCOND, and is almost always a slight 00184 * overestimate of the true error. 00185 * 00186 * BERR (output) REAL array, dimension (NRHS) 00187 * The componentwise relative backward error of each solution 00188 * vector X(j) (i.e., the smallest relative change in 00189 * any element of A or B that makes X(j) an exact solution). 00190 * 00191 * WORK (workspace) REAL array, dimension (3*N) 00192 * 00193 * IWORK (workspace) INTEGER array, dimension (N) 00194 * 00195 * INFO (output) INTEGER 00196 * = 0: successful exit 00197 * < 0: if INFO = -i, the i-th argument had an illegal value 00198 * > 0: if INFO = i, and i is 00199 * <= N: the leading minor of order i of A is 00200 * not positive definite, so the factorization 00201 * could not be completed, and the solution has not 00202 * been computed. RCOND = 0 is returned. 00203 * = N+1: U is nonsingular, but RCOND is less than machine 00204 * precision, meaning that the matrix is singular 00205 * to working precision. Nevertheless, the 00206 * solution and error bounds are computed because 00207 * there are a number of situations where the 00208 * computed solution can be more accurate than the 00209 * value of RCOND would suggest. 00210 * 00211 * ===================================================================== 00212 * 00213 * .. Parameters .. 00214 REAL ZERO, ONE 00215 PARAMETER ( ZERO = 0.0E+0, ONE = 1.0E+0 ) 00216 * .. 00217 * .. Local Scalars .. 00218 LOGICAL EQUIL, NOFACT, RCEQU 00219 INTEGER I, INFEQU, J 00220 REAL AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM 00221 * .. 00222 * .. External Functions .. 00223 LOGICAL LSAME 00224 REAL SLAMCH, SLANSY 00225 EXTERNAL LSAME, SLAMCH, SLANSY 00226 * .. 00227 * .. External Subroutines .. 00228 EXTERNAL SLACPY, SLAQSY, SPOCON, SPOEQU, SPORFS, SPOTRF, 00229 $ SPOTRS, XERBLA 00230 * .. 00231 * .. Intrinsic Functions .. 00232 INTRINSIC MAX, MIN 00233 * .. 00234 * .. Executable Statements .. 00235 * 00236 INFO = 0 00237 NOFACT = LSAME( FACT, 'N' ) 00238 EQUIL = LSAME( FACT, 'E' ) 00239 IF( NOFACT .OR. EQUIL ) THEN 00240 EQUED = 'N' 00241 RCEQU = .FALSE. 00242 ELSE 00243 RCEQU = LSAME( EQUED, 'Y' ) 00244 SMLNUM = SLAMCH( 'Safe minimum' ) 00245 BIGNUM = ONE / SMLNUM 00246 END IF 00247 * 00248 * Test the input parameters. 00249 * 00250 IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.LSAME( FACT, 'F' ) ) 00251 $ THEN 00252 INFO = -1 00253 ELSE IF( .NOT.LSAME( UPLO, 'U' ) .AND. .NOT.LSAME( UPLO, 'L' ) ) 00254 $ THEN 00255 INFO = -2 00256 ELSE IF( N.LT.0 ) THEN 00257 INFO = -3 00258 ELSE IF( NRHS.LT.0 ) THEN 00259 INFO = -4 00260 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN 00261 INFO = -6 00262 ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN 00263 INFO = -8 00264 ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT. 00265 $ ( RCEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN 00266 INFO = -9 00267 ELSE 00268 IF( RCEQU ) THEN 00269 SMIN = BIGNUM 00270 SMAX = ZERO 00271 DO 10 J = 1, N 00272 SMIN = MIN( SMIN, S( J ) ) 00273 SMAX = MAX( SMAX, S( J ) ) 00274 10 CONTINUE 00275 IF( SMIN.LE.ZERO ) THEN 00276 INFO = -10 00277 ELSE IF( N.GT.0 ) THEN 00278 SCOND = MAX( SMIN, SMLNUM ) / MIN( SMAX, BIGNUM ) 00279 ELSE 00280 SCOND = ONE 00281 END IF 00282 END IF 00283 IF( INFO.EQ.0 ) THEN 00284 IF( LDB.LT.MAX( 1, N ) ) THEN 00285 INFO = -12 00286 ELSE IF( LDX.LT.MAX( 1, N ) ) THEN 00287 INFO = -14 00288 END IF 00289 END IF 00290 END IF 00291 * 00292 IF( INFO.NE.0 ) THEN 00293 CALL XERBLA( 'SPOSVX', -INFO ) 00294 RETURN 00295 END IF 00296 * 00297 IF( EQUIL ) THEN 00298 * 00299 * Compute row and column scalings to equilibrate the matrix A. 00300 * 00301 CALL SPOEQU( N, A, LDA, S, SCOND, AMAX, INFEQU ) 00302 IF( INFEQU.EQ.0 ) THEN 00303 * 00304 * Equilibrate the matrix. 00305 * 00306 CALL SLAQSY( UPLO, N, A, LDA, S, SCOND, AMAX, EQUED ) 00307 RCEQU = LSAME( EQUED, 'Y' ) 00308 END IF 00309 END IF 00310 * 00311 * Scale the right hand side. 00312 * 00313 IF( RCEQU ) THEN 00314 DO 30 J = 1, NRHS 00315 DO 20 I = 1, N 00316 B( I, J ) = S( I )*B( I, J ) 00317 20 CONTINUE 00318 30 CONTINUE 00319 END IF 00320 * 00321 IF( NOFACT .OR. EQUIL ) THEN 00322 * 00323 * Compute the Cholesky factorization A = U**T *U or A = L*L**T. 00324 * 00325 CALL SLACPY( UPLO, N, N, A, LDA, AF, LDAF ) 00326 CALL SPOTRF( UPLO, N, AF, LDAF, INFO ) 00327 * 00328 * Return if INFO is non-zero. 00329 * 00330 IF( INFO.GT.0 )THEN 00331 RCOND = ZERO 00332 RETURN 00333 END IF 00334 END IF 00335 * 00336 * Compute the norm of the matrix A. 00337 * 00338 ANORM = SLANSY( '1', UPLO, N, A, LDA, WORK ) 00339 * 00340 * Compute the reciprocal of the condition number of A. 00341 * 00342 CALL SPOCON( UPLO, N, AF, LDAF, ANORM, RCOND, WORK, IWORK, INFO ) 00343 * 00344 * Compute the solution matrix X. 00345 * 00346 CALL SLACPY( 'Full', N, NRHS, B, LDB, X, LDX ) 00347 CALL SPOTRS( UPLO, N, NRHS, AF, LDAF, X, LDX, INFO ) 00348 * 00349 * Use iterative refinement to improve the computed solution and 00350 * compute error bounds and backward error estimates for it. 00351 * 00352 CALL SPORFS( UPLO, N, NRHS, A, LDA, AF, LDAF, B, LDB, X, LDX, 00353 $ FERR, BERR, WORK, IWORK, INFO ) 00354 * 00355 * Transform the solution matrix X to a solution of the original 00356 * system. 00357 * 00358 IF( RCEQU ) THEN 00359 DO 50 J = 1, NRHS 00360 DO 40 I = 1, N 00361 X( I, J ) = S( I )*X( I, J ) 00362 40 CONTINUE 00363 50 CONTINUE 00364 DO 60 J = 1, NRHS 00365 FERR( J ) = FERR( J ) / SCOND 00366 60 CONTINUE 00367 END IF 00368 * 00369 * Set INFO = N+1 if the matrix is singular to working precision. 00370 * 00371 IF( RCOND.LT.SLAMCH( 'Epsilon' ) ) 00372 $ INFO = N + 1 00373 * 00374 RETURN 00375 * 00376 * End of SPOSVX 00377 * 00378 END