![]() |
LAPACK
3.4.1
LAPACK: Linear Algebra PACKage
|
00001 *> \brief <b> CPPSVX computes the solution to system of linear equations A * X = B for OTHER matrices</b> 00002 * 00003 * =========== DOCUMENTATION =========== 00004 * 00005 * Online html documentation available at 00006 * http://www.netlib.org/lapack/explore-html/ 00007 * 00008 *> \htmlonly 00009 *> Download CPPSVX + dependencies 00010 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cppsvx.f"> 00011 *> [TGZ]</a> 00012 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cppsvx.f"> 00013 *> [ZIP]</a> 00014 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cppsvx.f"> 00015 *> [TXT]</a> 00016 *> \endhtmlonly 00017 * 00018 * Definition: 00019 * =========== 00020 * 00021 * SUBROUTINE CPPSVX( FACT, UPLO, N, NRHS, AP, AFP, EQUED, S, B, LDB, 00022 * X, LDX, RCOND, FERR, BERR, WORK, RWORK, INFO ) 00023 * 00024 * .. Scalar Arguments .. 00025 * CHARACTER EQUED, FACT, UPLO 00026 * INTEGER INFO, LDB, LDX, N, NRHS 00027 * REAL RCOND 00028 * .. 00029 * .. Array Arguments .. 00030 * REAL BERR( * ), FERR( * ), RWORK( * ), S( * ) 00031 * COMPLEX AFP( * ), AP( * ), B( LDB, * ), WORK( * ), 00032 * $ X( LDX, * ) 00033 * .. 00034 * 00035 * 00036 *> \par Purpose: 00037 * ============= 00038 *> 00039 *> \verbatim 00040 *> 00041 *> CPPSVX uses the Cholesky factorization A = U**H*U or A = L*L**H to 00042 *> compute the solution to a complex system of linear equations 00043 *> A * X = B, 00044 *> where A is an N-by-N Hermitian positive definite matrix stored in 00045 *> packed format and X and B are N-by-NRHS matrices. 00046 *> 00047 *> Error bounds on the solution and a condition estimate are also 00048 *> provided. 00049 *> \endverbatim 00050 * 00051 *> \par Description: 00052 * ================= 00053 *> 00054 *> \verbatim 00055 *> 00056 *> The following steps are performed: 00057 *> 00058 *> 1. If FACT = 'E', real scaling factors are computed to equilibrate 00059 *> the system: 00060 *> diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * B 00061 *> Whether or not the system will be equilibrated depends on the 00062 *> scaling of the matrix A, but if equilibration is used, A is 00063 *> overwritten by diag(S)*A*diag(S) and B by diag(S)*B. 00064 *> 00065 *> 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to 00066 *> factor the matrix A (after equilibration if FACT = 'E') as 00067 *> A = U**H * U , if UPLO = 'U', or 00068 *> A = L * L**H, if UPLO = 'L', 00069 *> where U is an upper triangular matrix, L is a lower triangular 00070 *> matrix, and **H indicates conjugate transpose. 00071 *> 00072 *> 3. If the leading i-by-i principal minor is not positive definite, 00073 *> then the routine returns with INFO = i. Otherwise, the factored 00074 *> form of A is used to estimate the condition number of the matrix 00075 *> A. If the reciprocal of the condition number is less than machine 00076 *> precision, INFO = N+1 is returned as a warning, but the routine 00077 *> still goes on to solve for X and compute error bounds as 00078 *> described below. 00079 *> 00080 *> 4. The system of equations is solved for X using the factored form 00081 *> of A. 00082 *> 00083 *> 5. Iterative refinement is applied to improve the computed solution 00084 *> matrix and calculate error bounds and backward error estimates 00085 *> for it. 00086 *> 00087 *> 6. If equilibration was used, the matrix X is premultiplied by 00088 *> diag(S) so that it solves the original system before 00089 *> equilibration. 00090 *> \endverbatim 00091 * 00092 * Arguments: 00093 * ========== 00094 * 00095 *> \param[in] FACT 00096 *> \verbatim 00097 *> FACT is CHARACTER*1 00098 *> Specifies whether or not the factored form of the matrix A is 00099 *> supplied on entry, and if not, whether the matrix A should be 00100 *> equilibrated before it is factored. 00101 *> = 'F': On entry, AFP contains the factored form of A. 00102 *> If EQUED = 'Y', the matrix A has been equilibrated 00103 *> with scaling factors given by S. AP and AFP will not 00104 *> be modified. 00105 *> = 'N': The matrix A will be copied to AFP and factored. 00106 *> = 'E': The matrix A will be equilibrated if necessary, then 00107 *> copied to AFP and factored. 00108 *> \endverbatim 00109 *> 00110 *> \param[in] UPLO 00111 *> \verbatim 00112 *> UPLO is CHARACTER*1 00113 *> = 'U': Upper triangle of A is stored; 00114 *> = 'L': Lower triangle of A is stored. 00115 *> \endverbatim 00116 *> 00117 *> \param[in] N 00118 *> \verbatim 00119 *> N is INTEGER 00120 *> The number of linear equations, i.e., the order of the 00121 *> matrix A. N >= 0. 00122 *> \endverbatim 00123 *> 00124 *> \param[in] NRHS 00125 *> \verbatim 00126 *> NRHS is INTEGER 00127 *> The number of right hand sides, i.e., the number of columns 00128 *> of the matrices B and X. NRHS >= 0. 00129 *> \endverbatim 00130 *> 00131 *> \param[in,out] AP 00132 *> \verbatim 00133 *> AP is COMPLEX array, dimension (N*(N+1)/2) 00134 *> On entry, the upper or lower triangle of the Hermitian matrix 00135 *> A, packed columnwise in a linear array, except if FACT = 'F' 00136 *> and EQUED = 'Y', then A must contain the equilibrated matrix 00137 *> diag(S)*A*diag(S). The j-th column of A is stored in the 00138 *> array AP as follows: 00139 *> if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j; 00140 *> if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n. 00141 *> See below for further details. A is not modified if 00142 *> FACT = 'F' or 'N', or if FACT = 'E' and EQUED = 'N' on exit. 00143 *> 00144 *> On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by 00145 *> diag(S)*A*diag(S). 00146 *> \endverbatim 00147 *> 00148 *> \param[in,out] AFP 00149 *> \verbatim 00150 *> AFP is COMPLEX array, dimension (N*(N+1)/2) 00151 *> If FACT = 'F', then AFP is an input argument and on entry 00152 *> contains the triangular factor U or L from the Cholesky 00153 *> factorization A = U**H*U or A = L*L**H, in the same storage 00154 *> format as A. If EQUED .ne. 'N', then AFP is the factored 00155 *> form of the equilibrated matrix A. 00156 *> 00157 *> If FACT = 'N', then AFP is an output argument and on exit 00158 *> returns the triangular factor U or L from the Cholesky 00159 *> factorization A = U**H * U or A = L * L**H of the original 00160 *> matrix A. 00161 *> 00162 *> If FACT = 'E', then AFP is an output argument and on exit 00163 *> returns the triangular factor U or L from the Cholesky 00164 *> factorization A = U**H*U or A = L*L**H of the equilibrated 00165 *> matrix A (see the description of AP for the form of the 00166 *> equilibrated matrix). 00167 *> \endverbatim 00168 *> 00169 *> \param[in,out] EQUED 00170 *> \verbatim 00171 *> EQUED is CHARACTER*1 00172 *> Specifies the form of equilibration that was done. 00173 *> = 'N': No equilibration (always true if FACT = 'N'). 00174 *> = 'Y': Equilibration was done, i.e., A has been replaced by 00175 *> diag(S) * A * diag(S). 00176 *> EQUED is an input argument if FACT = 'F'; otherwise, it is an 00177 *> output argument. 00178 *> \endverbatim 00179 *> 00180 *> \param[in,out] S 00181 *> \verbatim 00182 *> S is REAL array, dimension (N) 00183 *> The scale factors for A; not accessed if EQUED = 'N'. S is 00184 *> an input argument if FACT = 'F'; otherwise, S is an output 00185 *> argument. If FACT = 'F' and EQUED = 'Y', each element of S 00186 *> must be positive. 00187 *> \endverbatim 00188 *> 00189 *> \param[in,out] B 00190 *> \verbatim 00191 *> B is COMPLEX array, dimension (LDB,NRHS) 00192 *> On entry, the N-by-NRHS right hand side matrix B. 00193 *> On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y', 00194 *> B is overwritten by diag(S) * B. 00195 *> \endverbatim 00196 *> 00197 *> \param[in] LDB 00198 *> \verbatim 00199 *> LDB is INTEGER 00200 *> The leading dimension of the array B. LDB >= max(1,N). 00201 *> \endverbatim 00202 *> 00203 *> \param[out] X 00204 *> \verbatim 00205 *> X is COMPLEX array, dimension (LDX,NRHS) 00206 *> If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X to 00207 *> the original system of equations. Note that if EQUED = 'Y', 00208 *> A and B are modified on exit, and the solution to the 00209 *> equilibrated system is inv(diag(S))*X. 00210 *> \endverbatim 00211 *> 00212 *> \param[in] LDX 00213 *> \verbatim 00214 *> LDX is INTEGER 00215 *> The leading dimension of the array X. LDX >= max(1,N). 00216 *> \endverbatim 00217 *> 00218 *> \param[out] RCOND 00219 *> \verbatim 00220 *> RCOND is REAL 00221 *> The estimate of the reciprocal condition number of the matrix 00222 *> A after equilibration (if done). If RCOND is less than the 00223 *> machine precision (in particular, if RCOND = 0), the matrix 00224 *> is singular to working precision. This condition is 00225 *> indicated by a return code of INFO > 0. 00226 *> \endverbatim 00227 *> 00228 *> \param[out] FERR 00229 *> \verbatim 00230 *> FERR is REAL array, dimension (NRHS) 00231 *> The estimated forward error bound for each solution vector 00232 *> X(j) (the j-th column of the solution matrix X). 00233 *> If XTRUE is the true solution corresponding to X(j), FERR(j) 00234 *> is an estimated upper bound for the magnitude of the largest 00235 *> element in (X(j) - XTRUE) divided by the magnitude of the 00236 *> largest element in X(j). The estimate is as reliable as 00237 *> the estimate for RCOND, and is almost always a slight 00238 *> overestimate of the true error. 00239 *> \endverbatim 00240 *> 00241 *> \param[out] BERR 00242 *> \verbatim 00243 *> BERR is REAL array, dimension (NRHS) 00244 *> The componentwise relative backward error of each solution 00245 *> vector X(j) (i.e., the smallest relative change in 00246 *> any element of A or B that makes X(j) an exact solution). 00247 *> \endverbatim 00248 *> 00249 *> \param[out] WORK 00250 *> \verbatim 00251 *> WORK is COMPLEX array, dimension (2*N) 00252 *> \endverbatim 00253 *> 00254 *> \param[out] RWORK 00255 *> \verbatim 00256 *> RWORK is REAL array, dimension (N) 00257 *> \endverbatim 00258 *> 00259 *> \param[out] INFO 00260 *> \verbatim 00261 *> INFO is INTEGER 00262 *> = 0: successful exit 00263 *> < 0: if INFO = -i, the i-th argument had an illegal value 00264 *> > 0: if INFO = i, and i is 00265 *> <= N: the leading minor of order i of A is 00266 *> not positive definite, so the factorization 00267 *> could not be completed, and the solution has not 00268 *> been computed. RCOND = 0 is returned. 00269 *> = N+1: U is nonsingular, but RCOND is less than machine 00270 *> precision, meaning that the matrix is singular 00271 *> to working precision. Nevertheless, the 00272 *> solution and error bounds are computed because 00273 *> there are a number of situations where the 00274 *> computed solution can be more accurate than the 00275 *> value of RCOND would suggest. 00276 *> \endverbatim 00277 * 00278 * Authors: 00279 * ======== 00280 * 00281 *> \author Univ. of Tennessee 00282 *> \author Univ. of California Berkeley 00283 *> \author Univ. of Colorado Denver 00284 *> \author NAG Ltd. 00285 * 00286 *> \date April 2012 00287 * 00288 *> \ingroup complexOTHERsolve 00289 * 00290 *> \par Further Details: 00291 * ===================== 00292 *> 00293 *> \verbatim 00294 *> 00295 *> The packed storage scheme is illustrated by the following example 00296 *> when N = 4, UPLO = 'U': 00297 *> 00298 *> Two-dimensional storage of the Hermitian matrix A: 00299 *> 00300 *> a11 a12 a13 a14 00301 *> a22 a23 a24 00302 *> a33 a34 (aij = conjg(aji)) 00303 *> a44 00304 *> 00305 *> Packed storage of the upper triangle of A: 00306 *> 00307 *> AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ] 00308 *> \endverbatim 00309 *> 00310 * ===================================================================== 00311 SUBROUTINE CPPSVX( FACT, UPLO, N, NRHS, AP, AFP, EQUED, S, B, LDB, 00312 $ X, LDX, RCOND, FERR, BERR, WORK, RWORK, INFO ) 00313 * 00314 * -- LAPACK driver routine (version 3.4.1) -- 00315 * -- LAPACK is a software package provided by Univ. of Tennessee, -- 00316 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- 00317 * April 2012 00318 * 00319 * .. Scalar Arguments .. 00320 CHARACTER EQUED, FACT, UPLO 00321 INTEGER INFO, LDB, LDX, N, NRHS 00322 REAL RCOND 00323 * .. 00324 * .. Array Arguments .. 00325 REAL BERR( * ), FERR( * ), RWORK( * ), S( * ) 00326 COMPLEX AFP( * ), AP( * ), B( LDB, * ), WORK( * ), 00327 $ X( LDX, * ) 00328 * .. 00329 * 00330 * ===================================================================== 00331 * 00332 * .. Parameters .. 00333 REAL ZERO, ONE 00334 PARAMETER ( ZERO = 0.0E+0, ONE = 1.0E+0 ) 00335 * .. 00336 * .. Local Scalars .. 00337 LOGICAL EQUIL, NOFACT, RCEQU 00338 INTEGER I, INFEQU, J 00339 REAL AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM 00340 * .. 00341 * .. External Functions .. 00342 LOGICAL LSAME 00343 REAL CLANHP, SLAMCH 00344 EXTERNAL LSAME, CLANHP, SLAMCH 00345 * .. 00346 * .. External Subroutines .. 00347 EXTERNAL CCOPY, CLACPY, CLAQHP, CPPCON, CPPEQU, CPPRFS, 00348 $ CPPTRF, CPPTRS, XERBLA 00349 * .. 00350 * .. Intrinsic Functions .. 00351 INTRINSIC MAX, MIN 00352 * .. 00353 * .. Executable Statements .. 00354 * 00355 INFO = 0 00356 NOFACT = LSAME( FACT, 'N' ) 00357 EQUIL = LSAME( FACT, 'E' ) 00358 IF( NOFACT .OR. EQUIL ) THEN 00359 EQUED = 'N' 00360 RCEQU = .FALSE. 00361 ELSE 00362 RCEQU = LSAME( EQUED, 'Y' ) 00363 SMLNUM = SLAMCH( 'Safe minimum' ) 00364 BIGNUM = ONE / SMLNUM 00365 END IF 00366 * 00367 * Test the input parameters. 00368 * 00369 IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.LSAME( FACT, 'F' ) ) 00370 $ THEN 00371 INFO = -1 00372 ELSE IF( .NOT.LSAME( UPLO, 'U' ) .AND. .NOT.LSAME( UPLO, 'L' ) ) 00373 $ THEN 00374 INFO = -2 00375 ELSE IF( N.LT.0 ) THEN 00376 INFO = -3 00377 ELSE IF( NRHS.LT.0 ) THEN 00378 INFO = -4 00379 ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT. 00380 $ ( RCEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN 00381 INFO = -7 00382 ELSE 00383 IF( RCEQU ) THEN 00384 SMIN = BIGNUM 00385 SMAX = ZERO 00386 DO 10 J = 1, N 00387 SMIN = MIN( SMIN, S( J ) ) 00388 SMAX = MAX( SMAX, S( J ) ) 00389 10 CONTINUE 00390 IF( SMIN.LE.ZERO ) THEN 00391 INFO = -8 00392 ELSE IF( N.GT.0 ) THEN 00393 SCOND = MAX( SMIN, SMLNUM ) / MIN( SMAX, BIGNUM ) 00394 ELSE 00395 SCOND = ONE 00396 END IF 00397 END IF 00398 IF( INFO.EQ.0 ) THEN 00399 IF( LDB.LT.MAX( 1, N ) ) THEN 00400 INFO = -10 00401 ELSE IF( LDX.LT.MAX( 1, N ) ) THEN 00402 INFO = -12 00403 END IF 00404 END IF 00405 END IF 00406 * 00407 IF( INFO.NE.0 ) THEN 00408 CALL XERBLA( 'CPPSVX', -INFO ) 00409 RETURN 00410 END IF 00411 * 00412 IF( EQUIL ) THEN 00413 * 00414 * Compute row and column scalings to equilibrate the matrix A. 00415 * 00416 CALL CPPEQU( UPLO, N, AP, S, SCOND, AMAX, INFEQU ) 00417 IF( INFEQU.EQ.0 ) THEN 00418 * 00419 * Equilibrate the matrix. 00420 * 00421 CALL CLAQHP( UPLO, N, AP, S, SCOND, AMAX, EQUED ) 00422 RCEQU = LSAME( EQUED, 'Y' ) 00423 END IF 00424 END IF 00425 * 00426 * Scale the right-hand side. 00427 * 00428 IF( RCEQU ) THEN 00429 DO 30 J = 1, NRHS 00430 DO 20 I = 1, N 00431 B( I, J ) = S( I )*B( I, J ) 00432 20 CONTINUE 00433 30 CONTINUE 00434 END IF 00435 * 00436 IF( NOFACT .OR. EQUIL ) THEN 00437 * 00438 * Compute the Cholesky factorization A = U**H * U or A = L * L**H. 00439 * 00440 CALL CCOPY( N*( N+1 ) / 2, AP, 1, AFP, 1 ) 00441 CALL CPPTRF( UPLO, N, AFP, INFO ) 00442 * 00443 * Return if INFO is non-zero. 00444 * 00445 IF( INFO.GT.0 )THEN 00446 RCOND = ZERO 00447 RETURN 00448 END IF 00449 END IF 00450 * 00451 * Compute the norm of the matrix A. 00452 * 00453 ANORM = CLANHP( 'I', UPLO, N, AP, RWORK ) 00454 * 00455 * Compute the reciprocal of the condition number of A. 00456 * 00457 CALL CPPCON( UPLO, N, AFP, ANORM, RCOND, WORK, RWORK, INFO ) 00458 * 00459 * Compute the solution matrix X. 00460 * 00461 CALL CLACPY( 'Full', N, NRHS, B, LDB, X, LDX ) 00462 CALL CPPTRS( UPLO, N, NRHS, AFP, X, LDX, INFO ) 00463 * 00464 * Use iterative refinement to improve the computed solution and 00465 * compute error bounds and backward error estimates for it. 00466 * 00467 CALL CPPRFS( UPLO, N, NRHS, AP, AFP, B, LDB, X, LDX, FERR, BERR, 00468 $ WORK, RWORK, INFO ) 00469 * 00470 * Transform the solution matrix X to a solution of the original 00471 * system. 00472 * 00473 IF( RCEQU ) THEN 00474 DO 50 J = 1, NRHS 00475 DO 40 I = 1, N 00476 X( I, J ) = S( I )*X( I, J ) 00477 40 CONTINUE 00478 50 CONTINUE 00479 DO 60 J = 1, NRHS 00480 FERR( J ) = FERR( J ) / SCOND 00481 60 CONTINUE 00482 END IF 00483 * 00484 * Set INFO = N+1 if the matrix is singular to working precision. 00485 * 00486 IF( RCOND.LT.SLAMCH( 'Epsilon' ) ) 00487 $ INFO = N + 1 00488 * 00489 RETURN 00490 * 00491 * End of CPPSVX 00492 * 00493 END