LAPACK 3.3.1
Linear Algebra PACKage
|
00001 SUBROUTINE DPBSVX( FACT, UPLO, N, KD, NRHS, AB, LDAB, AFB, LDAFB, 00002 $ EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR, 00003 $ WORK, 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, KD, LDAB, LDAFB, LDB, LDX, N, NRHS 00013 DOUBLE PRECISION RCOND 00014 * .. 00015 * .. Array Arguments .. 00016 INTEGER IWORK( * ) 00017 DOUBLE PRECISION AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ), 00018 $ BERR( * ), FERR( * ), S( * ), WORK( * ), 00019 $ X( LDX, * ) 00020 * .. 00021 * 00022 * Purpose 00023 * ======= 00024 * 00025 * DPBSVX 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 band matrix and X 00029 * and B 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 band matrix, and L is a lower 00051 * triangular band 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, AFB contains the factored form of A. 00080 * If EQUED = 'Y', the matrix A has been equilibrated 00081 * with scaling factors given by S. AB and AFB will not 00082 * be modified. 00083 * = 'N': The matrix A will be copied to AFB and factored. 00084 * = 'E': The matrix A will be equilibrated if necessary, then 00085 * copied to AFB 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 * KD (input) INTEGER 00096 * The number of superdiagonals of the matrix A if UPLO = 'U', 00097 * or the number of subdiagonals if UPLO = 'L'. KD >= 0. 00098 * 00099 * NRHS (input) INTEGER 00100 * The number of right-hand sides, i.e., the number of columns 00101 * of the matrices B and X. NRHS >= 0. 00102 * 00103 * AB (input/output) DOUBLE PRECISION array, dimension (LDAB,N) 00104 * On entry, the upper or lower triangle of the symmetric band 00105 * matrix A, stored in the first KD+1 rows of the array, except 00106 * if FACT = 'F' and EQUED = 'Y', then A must contain the 00107 * equilibrated matrix diag(S)*A*diag(S). The j-th column of A 00108 * is stored in the j-th column of the array AB as follows: 00109 * if UPLO = 'U', AB(KD+1+i-j,j) = A(i,j) for max(1,j-KD)<=i<=j; 00110 * if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(N,j+KD). 00111 * See below for further details. 00112 * 00113 * On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by 00114 * diag(S)*A*diag(S). 00115 * 00116 * LDAB (input) INTEGER 00117 * The leading dimension of the array A. LDAB >= KD+1. 00118 * 00119 * AFB (input or output) DOUBLE PRECISION array, dimension (LDAFB,N) 00120 * If FACT = 'F', then AFB is an input argument and on entry 00121 * contains the triangular factor U or L from the Cholesky 00122 * factorization A = U**T*U or A = L*L**T of the band matrix 00123 * A, in the same storage format as A (see AB). If EQUED = 'Y', 00124 * then AFB is the factored form of the equilibrated matrix A. 00125 * 00126 * If FACT = 'N', then AFB is an output argument and on exit 00127 * returns the triangular factor U or L from the Cholesky 00128 * factorization A = U**T*U or A = L*L**T. 00129 * 00130 * If FACT = 'E', then AFB is an output argument and on exit 00131 * returns the triangular factor U or L from the Cholesky 00132 * factorization A = U**T*U or A = L*L**T of the equilibrated 00133 * matrix A (see the description of A for the form of the 00134 * equilibrated matrix). 00135 * 00136 * LDAFB (input) INTEGER 00137 * The leading dimension of the array AFB. LDAFB >= KD+1. 00138 * 00139 * EQUED (input or output) CHARACTER*1 00140 * Specifies the form of equilibration that was done. 00141 * = 'N': No equilibration (always true if FACT = 'N'). 00142 * = 'Y': Equilibration was done, i.e., A has been replaced by 00143 * diag(S) * A * diag(S). 00144 * EQUED is an input argument if FACT = 'F'; otherwise, it is an 00145 * output argument. 00146 * 00147 * S (input or output) DOUBLE PRECISION array, dimension (N) 00148 * The scale factors for A; not accessed if EQUED = 'N'. S is 00149 * an input argument if FACT = 'F'; otherwise, S is an output 00150 * argument. If FACT = 'F' and EQUED = 'Y', each element of S 00151 * must be positive. 00152 * 00153 * B (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS) 00154 * On entry, the N-by-NRHS right hand side matrix B. 00155 * On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y', 00156 * B is overwritten by diag(S) * B. 00157 * 00158 * LDB (input) INTEGER 00159 * The leading dimension of the array B. LDB >= max(1,N). 00160 * 00161 * X (output) DOUBLE PRECISION array, dimension (LDX,NRHS) 00162 * If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X to 00163 * the original system of equations. Note that if EQUED = 'Y', 00164 * A and B are modified on exit, and the solution to the 00165 * equilibrated system is inv(diag(S))*X. 00166 * 00167 * LDX (input) INTEGER 00168 * The leading dimension of the array X. LDX >= max(1,N). 00169 * 00170 * RCOND (output) DOUBLE PRECISION 00171 * The estimate of the reciprocal condition number of the matrix 00172 * A after equilibration (if done). If RCOND is less than the 00173 * machine precision (in particular, if RCOND = 0), the matrix 00174 * is singular to working precision. This condition is 00175 * indicated by a return code of INFO > 0. 00176 * 00177 * FERR (output) DOUBLE PRECISION array, dimension (NRHS) 00178 * The estimated forward error bound for each solution vector 00179 * X(j) (the j-th column of the solution matrix X). 00180 * If XTRUE is the true solution corresponding to X(j), FERR(j) 00181 * is an estimated upper bound for the magnitude of the largest 00182 * element in (X(j) - XTRUE) divided by the magnitude of the 00183 * largest element in X(j). The estimate is as reliable as 00184 * the estimate for RCOND, and is almost always a slight 00185 * overestimate of the true error. 00186 * 00187 * BERR (output) DOUBLE PRECISION array, dimension (NRHS) 00188 * The componentwise relative backward error of each solution 00189 * vector X(j) (i.e., the smallest relative change in 00190 * any element of A or B that makes X(j) an exact solution). 00191 * 00192 * WORK (workspace) DOUBLE PRECISION array, dimension (3*N) 00193 * 00194 * IWORK (workspace) INTEGER array, dimension (N) 00195 * 00196 * INFO (output) INTEGER 00197 * = 0: successful exit 00198 * < 0: if INFO = -i, the i-th argument had an illegal value 00199 * > 0: if INFO = i, and i is 00200 * <= N: the leading minor of order i of A is 00201 * not positive definite, so the factorization 00202 * could not be completed, and the solution has not 00203 * been computed. RCOND = 0 is returned. 00204 * = N+1: U is nonsingular, but RCOND is less than machine 00205 * precision, meaning that the matrix is singular 00206 * to working precision. Nevertheless, the 00207 * solution and error bounds are computed because 00208 * there are a number of situations where the 00209 * computed solution can be more accurate than the 00210 * value of RCOND would suggest. 00211 * 00212 * Further Details 00213 * =============== 00214 * 00215 * The band storage scheme is illustrated by the following example, when 00216 * N = 6, KD = 2, and UPLO = 'U': 00217 * 00218 * Two-dimensional storage of the symmetric matrix A: 00219 * 00220 * a11 a12 a13 00221 * a22 a23 a24 00222 * a33 a34 a35 00223 * a44 a45 a46 00224 * a55 a56 00225 * (aij=conjg(aji)) a66 00226 * 00227 * Band storage of the upper triangle of A: 00228 * 00229 * * * a13 a24 a35 a46 00230 * * a12 a23 a34 a45 a56 00231 * a11 a22 a33 a44 a55 a66 00232 * 00233 * Similarly, if UPLO = 'L' the format of A is as follows: 00234 * 00235 * a11 a22 a33 a44 a55 a66 00236 * a21 a32 a43 a54 a65 * 00237 * a31 a42 a53 a64 * * 00238 * 00239 * Array elements marked * are not used by the routine. 00240 * 00241 * ===================================================================== 00242 * 00243 * .. Parameters .. 00244 DOUBLE PRECISION ZERO, ONE 00245 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 ) 00246 * .. 00247 * .. Local Scalars .. 00248 LOGICAL EQUIL, NOFACT, RCEQU, UPPER 00249 INTEGER I, INFEQU, J, J1, J2 00250 DOUBLE PRECISION AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM 00251 * .. 00252 * .. External Functions .. 00253 LOGICAL LSAME 00254 DOUBLE PRECISION DLAMCH, DLANSB 00255 EXTERNAL LSAME, DLAMCH, DLANSB 00256 * .. 00257 * .. External Subroutines .. 00258 EXTERNAL DCOPY, DLACPY, DLAQSB, DPBCON, DPBEQU, DPBRFS, 00259 $ DPBTRF, DPBTRS, XERBLA 00260 * .. 00261 * .. Intrinsic Functions .. 00262 INTRINSIC MAX, MIN 00263 * .. 00264 * .. Executable Statements .. 00265 * 00266 INFO = 0 00267 NOFACT = LSAME( FACT, 'N' ) 00268 EQUIL = LSAME( FACT, 'E' ) 00269 UPPER = LSAME( UPLO, 'U' ) 00270 IF( NOFACT .OR. EQUIL ) THEN 00271 EQUED = 'N' 00272 RCEQU = .FALSE. 00273 ELSE 00274 RCEQU = LSAME( EQUED, 'Y' ) 00275 SMLNUM = DLAMCH( 'Safe minimum' ) 00276 BIGNUM = ONE / SMLNUM 00277 END IF 00278 * 00279 * Test the input parameters. 00280 * 00281 IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.LSAME( FACT, 'F' ) ) 00282 $ THEN 00283 INFO = -1 00284 ELSE IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN 00285 INFO = -2 00286 ELSE IF( N.LT.0 ) THEN 00287 INFO = -3 00288 ELSE IF( KD.LT.0 ) THEN 00289 INFO = -4 00290 ELSE IF( NRHS.LT.0 ) THEN 00291 INFO = -5 00292 ELSE IF( LDAB.LT.KD+1 ) THEN 00293 INFO = -7 00294 ELSE IF( LDAFB.LT.KD+1 ) THEN 00295 INFO = -9 00296 ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT. 00297 $ ( RCEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN 00298 INFO = -10 00299 ELSE 00300 IF( RCEQU ) THEN 00301 SMIN = BIGNUM 00302 SMAX = ZERO 00303 DO 10 J = 1, N 00304 SMIN = MIN( SMIN, S( J ) ) 00305 SMAX = MAX( SMAX, S( J ) ) 00306 10 CONTINUE 00307 IF( SMIN.LE.ZERO ) THEN 00308 INFO = -11 00309 ELSE IF( N.GT.0 ) THEN 00310 SCOND = MAX( SMIN, SMLNUM ) / MIN( SMAX, BIGNUM ) 00311 ELSE 00312 SCOND = ONE 00313 END IF 00314 END IF 00315 IF( INFO.EQ.0 ) THEN 00316 IF( LDB.LT.MAX( 1, N ) ) THEN 00317 INFO = -13 00318 ELSE IF( LDX.LT.MAX( 1, N ) ) THEN 00319 INFO = -15 00320 END IF 00321 END IF 00322 END IF 00323 * 00324 IF( INFO.NE.0 ) THEN 00325 CALL XERBLA( 'DPBSVX', -INFO ) 00326 RETURN 00327 END IF 00328 * 00329 IF( EQUIL ) THEN 00330 * 00331 * Compute row and column scalings to equilibrate the matrix A. 00332 * 00333 CALL DPBEQU( UPLO, N, KD, AB, LDAB, S, SCOND, AMAX, INFEQU ) 00334 IF( INFEQU.EQ.0 ) THEN 00335 * 00336 * Equilibrate the matrix. 00337 * 00338 CALL DLAQSB( UPLO, N, KD, AB, LDAB, S, SCOND, AMAX, EQUED ) 00339 RCEQU = LSAME( EQUED, 'Y' ) 00340 END IF 00341 END IF 00342 * 00343 * Scale the right-hand side. 00344 * 00345 IF( RCEQU ) THEN 00346 DO 30 J = 1, NRHS 00347 DO 20 I = 1, N 00348 B( I, J ) = S( I )*B( I, J ) 00349 20 CONTINUE 00350 30 CONTINUE 00351 END IF 00352 * 00353 IF( NOFACT .OR. EQUIL ) THEN 00354 * 00355 * Compute the Cholesky factorization A = U**T *U or A = L*L**T. 00356 * 00357 IF( UPPER ) THEN 00358 DO 40 J = 1, N 00359 J1 = MAX( J-KD, 1 ) 00360 CALL DCOPY( J-J1+1, AB( KD+1-J+J1, J ), 1, 00361 $ AFB( KD+1-J+J1, J ), 1 ) 00362 40 CONTINUE 00363 ELSE 00364 DO 50 J = 1, N 00365 J2 = MIN( J+KD, N ) 00366 CALL DCOPY( J2-J+1, AB( 1, J ), 1, AFB( 1, J ), 1 ) 00367 50 CONTINUE 00368 END IF 00369 * 00370 CALL DPBTRF( UPLO, N, KD, AFB, LDAFB, INFO ) 00371 * 00372 * Return if INFO is non-zero. 00373 * 00374 IF( INFO.GT.0 )THEN 00375 RCOND = ZERO 00376 RETURN 00377 END IF 00378 END IF 00379 * 00380 * Compute the norm of the matrix A. 00381 * 00382 ANORM = DLANSB( '1', UPLO, N, KD, AB, LDAB, WORK ) 00383 * 00384 * Compute the reciprocal of the condition number of A. 00385 * 00386 CALL DPBCON( UPLO, N, KD, AFB, LDAFB, ANORM, RCOND, WORK, IWORK, 00387 $ INFO ) 00388 * 00389 * Compute the solution matrix X. 00390 * 00391 CALL DLACPY( 'Full', N, NRHS, B, LDB, X, LDX ) 00392 CALL DPBTRS( UPLO, N, KD, NRHS, AFB, LDAFB, X, LDX, INFO ) 00393 * 00394 * Use iterative refinement to improve the computed solution and 00395 * compute error bounds and backward error estimates for it. 00396 * 00397 CALL DPBRFS( UPLO, N, KD, NRHS, AB, LDAB, AFB, LDAFB, B, LDB, X, 00398 $ LDX, FERR, BERR, WORK, IWORK, INFO ) 00399 * 00400 * Transform the solution matrix X to a solution of the original 00401 * system. 00402 * 00403 IF( RCEQU ) THEN 00404 DO 70 J = 1, NRHS 00405 DO 60 I = 1, N 00406 X( I, J ) = S( I )*X( I, J ) 00407 60 CONTINUE 00408 70 CONTINUE 00409 DO 80 J = 1, NRHS 00410 FERR( J ) = FERR( J ) / SCOND 00411 80 CONTINUE 00412 END IF 00413 * 00414 * Set INFO = N+1 if the matrix is singular to working precision. 00415 * 00416 IF( RCOND.LT.DLAMCH( 'Epsilon' ) ) 00417 $ INFO = N + 1 00418 * 00419 RETURN 00420 * 00421 * End of DPBSVX 00422 * 00423 END