1 *> \brief <b> ZPOSVX 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 ZPOSVX + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zposvx.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zposvx.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zposvx.f">
21 * SUBROUTINE ZPOSVX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED,
22 * S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK,
25 * .. Scalar Arguments ..
26 * CHARACTER EQUED, FACT, UPLO
27 * INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS
28 * DOUBLE PRECISION RCOND
30 * .. Array Arguments ..
31 * DOUBLE PRECISION BERR( * ), FERR( * ), RWORK( * ), S( * )
32 * COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
33 * $ WORK( * ), X( LDX, * )
42 *> ZPOSVX uses the Cholesky factorization A = U**H*U or A = L*L**H to
43 *> compute the solution to a complex system of linear equations
45 *> where A is an N-by-N Hermitian positive definite matrix and X and B
46 *> are N-by-NRHS matrices.
48 *> Error bounds on the solution and a condition estimate are also
57 *> The following steps are performed:
59 *> 1. If FACT = 'E', real scaling factors are computed to equilibrate
61 *> diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * B
62 *> Whether or not the system will be equilibrated depends on the
63 *> scaling of the matrix A, but if equilibration is used, A is
64 *> overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
66 *> 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
67 *> factor the matrix A (after equilibration if FACT = 'E') as
68 *> A = U**H* U, if UPLO = 'U', or
69 *> A = L * L**H, if UPLO = 'L',
70 *> where U is an upper triangular matrix and L is a lower triangular
73 *> 3. If the leading i-by-i principal minor is not positive definite,
74 *> then the routine returns with INFO = i. Otherwise, the factored
75 *> form of A is used to estimate the condition number of the matrix
76 *> A. If the reciprocal of the condition number is less than machine
77 *> precision, INFO = N+1 is returned as a warning, but the routine
78 *> still goes on to solve for X and compute error bounds as
81 *> 4. The system of equations is solved for X using the factored form
84 *> 5. Iterative refinement is applied to improve the computed solution
85 *> matrix and calculate error bounds and backward error estimates
88 *> 6. If equilibration was used, the matrix X is premultiplied by
89 *> diag(S) so that it solves the original system before
98 *> FACT is CHARACTER*1
99 *> Specifies whether or not the factored form of the matrix A is
100 *> supplied on entry, and if not, whether the matrix A should be
101 *> equilibrated before it is factored.
102 *> = 'F': On entry, AF contains the factored form of A.
103 *> If EQUED = 'Y', the matrix A has been equilibrated
104 *> with scaling factors given by S. A and AF will not
106 *> = 'N': The matrix A will be copied to AF and factored.
107 *> = 'E': The matrix A will be equilibrated if necessary, then
108 *> copied to AF and factored.
113 *> UPLO is CHARACTER*1
114 *> = 'U': Upper triangle of A is stored;
115 *> = 'L': Lower triangle of A is stored.
121 *> The number of linear equations, i.e., the order of the
128 *> The number of right hand sides, i.e., the number of columns
129 *> of the matrices B and X. NRHS >= 0.
134 *> A is COMPLEX*16 array, dimension (LDA,N)
135 *> On entry, the Hermitian matrix A, except if FACT = 'F' and
136 *> EQUED = 'Y', then A must contain the equilibrated matrix
137 *> diag(S)*A*diag(S). If UPLO = 'U', the leading
138 *> N-by-N upper triangular part of A contains the upper
139 *> triangular part of the matrix A, and the strictly lower
140 *> triangular part of A is not referenced. If UPLO = 'L', the
141 *> leading N-by-N lower triangular part of A contains the lower
142 *> triangular part of the matrix A, and the strictly upper
143 *> triangular part of A is not referenced. A is not modified if
144 *> FACT = 'F' or 'N', or if FACT = 'E' and EQUED = 'N' on exit.
146 *> On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
147 *> diag(S)*A*diag(S).
153 *> The leading dimension of the array A. LDA >= max(1,N).
158 *> AF is COMPLEX*16 array, dimension (LDAF,N)
159 *> If FACT = 'F', then AF is an input argument and on entry
160 *> contains the triangular factor U or L from the Cholesky
161 *> factorization A = U**H *U or A = L*L**H, in the same storage
162 *> format as A. If EQUED .ne. 'N', then AF is the factored form
163 *> of the equilibrated matrix diag(S)*A*diag(S).
165 *> If FACT = 'N', then AF is an output argument and on exit
166 *> returns the triangular factor U or L from the Cholesky
167 *> factorization A = U**H *U or A = L*L**H of the original
170 *> If FACT = 'E', then AF is an output argument and on exit
171 *> returns the triangular factor U or L from the Cholesky
172 *> factorization A = U**H *U or A = L*L**H of the equilibrated
173 *> matrix A (see the description of A for the form of the
174 *> equilibrated matrix).
180 *> The leading dimension of the array AF. LDAF >= max(1,N).
183 *> \param[in,out] EQUED
185 *> EQUED is CHARACTER*1
186 *> Specifies the form of equilibration that was done.
187 *> = 'N': No equilibration (always true if FACT = 'N').
188 *> = 'Y': Equilibration was done, i.e., A has been replaced by
189 *> diag(S) * A * diag(S).
190 *> EQUED is an input argument if FACT = 'F'; otherwise, it is an
196 *> S is DOUBLE PRECISION array, dimension (N)
197 *> The scale factors for A; not accessed if EQUED = 'N'. S is
198 *> an input argument if FACT = 'F'; otherwise, S is an output
199 *> argument. If FACT = 'F' and EQUED = 'Y', each element of S
205 *> B is COMPLEX*16 array, dimension (LDB,NRHS)
206 *> On entry, the N-by-NRHS righthand side matrix B.
207 *> On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y',
208 *> B is overwritten by diag(S) * B.
214 *> The leading dimension of the array B. LDB >= max(1,N).
219 *> X is COMPLEX*16 array, dimension (LDX,NRHS)
220 *> If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X to
221 *> the original system of equations. Note that if EQUED = 'Y',
222 *> A and B are modified on exit, and the solution to the
223 *> equilibrated system is inv(diag(S))*X.
229 *> The leading dimension of the array X. LDX >= max(1,N).
234 *> RCOND is DOUBLE PRECISION
235 *> The estimate of the reciprocal condition number of the matrix
236 *> A after equilibration (if done). If RCOND is less than the
237 *> machine precision (in particular, if RCOND = 0), the matrix
238 *> is singular to working precision. This condition is
239 *> indicated by a return code of INFO > 0.
244 *> FERR is DOUBLE PRECISION array, dimension (NRHS)
245 *> The estimated forward error bound for each solution vector
246 *> X(j) (the j-th column of the solution matrix X).
247 *> If XTRUE is the true solution corresponding to X(j), FERR(j)
248 *> is an estimated upper bound for the magnitude of the largest
249 *> element in (X(j) - XTRUE) divided by the magnitude of the
250 *> largest element in X(j). The estimate is as reliable as
251 *> the estimate for RCOND, and is almost always a slight
252 *> overestimate of the true error.
257 *> BERR is DOUBLE PRECISION array, dimension (NRHS)
258 *> The componentwise relative backward error of each solution
259 *> vector X(j) (i.e., the smallest relative change in
260 *> any element of A or B that makes X(j) an exact solution).
265 *> WORK is COMPLEX*16 array, dimension (2*N)
270 *> RWORK is DOUBLE PRECISION array, dimension (N)
276 *> = 0: successful exit
277 *> < 0: if INFO = -i, the i-th argument had an illegal value
278 *> > 0: if INFO = i, and i is
279 *> <= N: the leading minor of order i of A is
280 *> not positive definite, so the factorization
281 *> could not be completed, and the solution has not
282 *> been computed. RCOND = 0 is returned.
283 *> = N+1: U is nonsingular, but RCOND is less than machine
284 *> precision, meaning that the matrix is singular
285 *> to working precision. Nevertheless, the
286 *> solution and error bounds are computed because
287 *> there are a number of situations where the
288 *> computed solution can be more accurate than the
289 *> value of RCOND would suggest.
295 *> \author Univ. of Tennessee
296 *> \author Univ. of California Berkeley
297 *> \author Univ. of Colorado Denver
302 *> \ingroup complex16POsolve
304 * =====================================================================
305 SUBROUTINE ZPOSVX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED,
306 $ S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK,
309 * -- LAPACK driver routine (version 3.4.1) --
310 * -- LAPACK is a software package provided by Univ. of Tennessee, --
311 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
314 * .. Scalar Arguments ..
315 CHARACTER EQUED, FACT, UPLO
316 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS
317 DOUBLE PRECISION RCOND
319 * .. Array Arguments ..
320 DOUBLE PRECISION BERR( * ), FERR( * ), RWORK( * ), S( * )
321 COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
322 $ WORK( * ), X( LDX, * )
325 * =====================================================================
328 DOUBLE PRECISION ZERO, ONE
329 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
331 * .. Local Scalars ..
332 LOGICAL EQUIL, NOFACT, RCEQU
334 DOUBLE PRECISION AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM
336 * .. External Functions ..
338 DOUBLE PRECISION DLAMCH, ZLANHE
339 EXTERNAL LSAME, DLAMCH, ZLANHE
341 * .. External Subroutines ..
342 EXTERNAL XERBLA, ZLACPY, ZLAQHE, ZPOCON, ZPOEQU, ZPORFS,
345 * .. Intrinsic Functions ..
348 * .. Executable Statements ..
351 NOFACT = LSAME( FACT, 'N' )
352 EQUIL = LSAME( FACT, 'E' )
353 IF( NOFACT .OR. EQUIL ) THEN
357 RCEQU = LSAME( EQUED, 'Y' )
358 SMLNUM = DLAMCH( 'Safe minimum' )
359 BIGNUM = ONE / SMLNUM
362 * Test the input parameters.
364 IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.LSAME( FACT, 'F' ) )
367 ELSE IF( .NOT.LSAME( UPLO, 'U' ) .AND. .NOT.LSAME( UPLO, 'L' ) )
370 ELSE IF( N.LT.0 ) THEN
372 ELSE IF( NRHS.LT.0 ) THEN
374 ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
376 ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN
378 ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
379 $ ( RCEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
386 SMIN = MIN( SMIN, S( J ) )
387 SMAX = MAX( SMAX, S( J ) )
389 IF( SMIN.LE.ZERO ) THEN
391 ELSE IF( N.GT.0 ) THEN
392 SCOND = MAX( SMIN, SMLNUM ) / MIN( SMAX, BIGNUM )
398 IF( LDB.LT.MAX( 1, N ) ) THEN
400 ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
407 CALL XERBLA( 'ZPOSVX', -INFO )
413 * Compute row and column scalings to equilibrate the matrix A.
415 CALL ZPOEQU( N, A, LDA, S, SCOND, AMAX, INFEQU )
416 IF( INFEQU.EQ.0 ) THEN
418 * Equilibrate the matrix.
420 CALL ZLAQHE( UPLO, N, A, LDA, S, SCOND, AMAX, EQUED )
421 RCEQU = LSAME( EQUED, 'Y' )
425 * Scale the right hand side.
430 B( I, J ) = S( I )*B( I, J )
435 IF( NOFACT .OR. EQUIL ) THEN
437 * Compute the Cholesky factorization A = U**H *U or A = L*L**H.
439 CALL ZLACPY( UPLO, N, N, A, LDA, AF, LDAF )
440 CALL ZPOTRF( UPLO, N, AF, LDAF, INFO )
442 * Return if INFO is non-zero.
450 * Compute the norm of the matrix A.
452 ANORM = ZLANHE( '1', UPLO, N, A, LDA, RWORK )
454 * Compute the reciprocal of the condition number of A.
456 CALL ZPOCON( UPLO, N, AF, LDAF, ANORM, RCOND, WORK, RWORK, INFO )
458 * Compute the solution matrix X.
460 CALL ZLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
461 CALL ZPOTRS( UPLO, N, NRHS, AF, LDAF, X, LDX, INFO )
463 * Use iterative refinement to improve the computed solution and
464 * compute error bounds and backward error estimates for it.
466 CALL ZPORFS( UPLO, N, NRHS, A, LDA, AF, LDAF, B, LDB, X, LDX,
467 $ FERR, BERR, WORK, RWORK, INFO )
469 * Transform the solution matrix X to a solution of the original
475 X( I, J ) = S( I )*X( I, J )
479 FERR( J ) = FERR( J ) / SCOND
483 * Set INFO = N+1 if the matrix is singular to working precision.
485 IF( RCOND.LT.DLAMCH( 'Epsilon' ) )