LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zposvx.f
Go to the documentation of this file.
1*> \brief <b> ZPOSVX computes the solution to system of linear equations A * X = B for PO matrices</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download ZPOSVX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zposvx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zposvx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zposvx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE ZPOSVX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED,
20* S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK,
21* RWORK, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER EQUED, FACT, UPLO
25* INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS
26* DOUBLE PRECISION RCOND
27* ..
28* .. Array Arguments ..
29* DOUBLE PRECISION BERR( * ), FERR( * ), RWORK( * ), S( * )
30* COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
31* $ WORK( * ), X( LDX, * )
32* ..
33*
34*
35*> \par Purpose:
36* =============
37*>
38*> \verbatim
39*>
40*> ZPOSVX uses the Cholesky factorization A = U**H*U or A = L*L**H to
41*> compute the solution to a complex system of linear equations
42*> A * X = B,
43*> where A is an N-by-N Hermitian positive definite matrix and X and B
44*> are N-by-NRHS matrices.
45*>
46*> Error bounds on the solution and a condition estimate are also
47*> provided.
48*> \endverbatim
49*
50*> \par Description:
51* =================
52*>
53*> \verbatim
54*>
55*> The following steps are performed:
56*>
57*> 1. If FACT = 'E', real scaling factors are computed to equilibrate
58*> the system:
59*> diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * B
60*> Whether or not the system will be equilibrated depends on the
61*> scaling of the matrix A, but if equilibration is used, A is
62*> overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
63*>
64*> 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
65*> factor the matrix A (after equilibration if FACT = 'E') as
66*> A = U**H* U, if UPLO = 'U', or
67*> A = L * L**H, if UPLO = 'L',
68*> where U is an upper triangular matrix and L is a lower triangular
69*> matrix.
70*>
71*> 3. If the leading principal minor of order i is not positive,
72*> then the routine returns with INFO = i. Otherwise, the factored
73*> form of A is used to estimate the condition number of the matrix
74*> A. If the reciprocal of the condition number is less than machine
75*> precision, INFO = N+1 is returned as a warning, but the routine
76*> still goes on to solve for X and compute error bounds as
77*> described below.
78*>
79*> 4. The system of equations is solved for X using the factored form
80*> of A.
81*>
82*> 5. Iterative refinement is applied to improve the computed solution
83*> matrix and calculate error bounds and backward error estimates
84*> for it.
85*>
86*> 6. If equilibration was used, the matrix X is premultiplied by
87*> diag(S) so that it solves the original system before
88*> equilibration.
89*> \endverbatim
90*
91* Arguments:
92* ==========
93*
94*> \param[in] FACT
95*> \verbatim
96*> FACT is CHARACTER*1
97*> Specifies whether or not the factored form of the matrix A is
98*> supplied on entry, and if not, whether the matrix A should be
99*> equilibrated before it is factored.
100*> = 'F': On entry, AF contains the factored form of A.
101*> If EQUED = 'Y', the matrix A has been equilibrated
102*> with scaling factors given by S. A and AF will not
103*> be modified.
104*> = 'N': The matrix A will be copied to AF and factored.
105*> = 'E': The matrix A will be equilibrated if necessary, then
106*> copied to AF and factored.
107*> \endverbatim
108*>
109*> \param[in] UPLO
110*> \verbatim
111*> UPLO is CHARACTER*1
112*> = 'U': Upper triangle of A is stored;
113*> = 'L': Lower triangle of A is stored.
114*> \endverbatim
115*>
116*> \param[in] N
117*> \verbatim
118*> N is INTEGER
119*> The number of linear equations, i.e., the order of the
120*> matrix A. N >= 0.
121*> \endverbatim
122*>
123*> \param[in] NRHS
124*> \verbatim
125*> NRHS is INTEGER
126*> The number of right hand sides, i.e., the number of columns
127*> of the matrices B and X. NRHS >= 0.
128*> \endverbatim
129*>
130*> \param[in,out] A
131*> \verbatim
132*> A is COMPLEX*16 array, dimension (LDA,N)
133*> On entry, the Hermitian matrix A, except if FACT = 'F' and
134*> EQUED = 'Y', then A must contain the equilibrated matrix
135*> diag(S)*A*diag(S). If UPLO = 'U', the leading
136*> N-by-N upper triangular part of A contains the upper
137*> triangular part of the matrix A, and the strictly lower
138*> triangular part of A is not referenced. If UPLO = 'L', the
139*> leading N-by-N lower triangular part of A contains the lower
140*> triangular part of the matrix A, and the strictly upper
141*> triangular part of A is not referenced. A is not modified if
142*> FACT = 'F' or 'N', or if FACT = 'E' and EQUED = 'N' on exit.
143*>
144*> On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
145*> diag(S)*A*diag(S).
146*> \endverbatim
147*>
148*> \param[in] LDA
149*> \verbatim
150*> LDA is INTEGER
151*> The leading dimension of the array A. LDA >= max(1,N).
152*> \endverbatim
153*>
154*> \param[in,out] AF
155*> \verbatim
156*> AF is COMPLEX*16 array, dimension (LDAF,N)
157*> If FACT = 'F', then AF is an input argument and on entry
158*> contains the triangular factor U or L from the Cholesky
159*> factorization A = U**H *U or A = L*L**H, in the same storage
160*> format as A. If EQUED .ne. 'N', then AF is the factored form
161*> of the equilibrated matrix diag(S)*A*diag(S).
162*>
163*> If FACT = 'N', then AF is an output argument and on exit
164*> returns the triangular factor U or L from the Cholesky
165*> factorization A = U**H *U or A = L*L**H of the original
166*> matrix A.
167*>
168*> If FACT = 'E', then AF is an output argument and on exit
169*> returns the triangular factor U or L from the Cholesky
170*> factorization A = U**H *U or A = L*L**H of the equilibrated
171*> matrix A (see the description of A for the form of the
172*> equilibrated matrix).
173*> \endverbatim
174*>
175*> \param[in] LDAF
176*> \verbatim
177*> LDAF is INTEGER
178*> The leading dimension of the array AF. LDAF >= max(1,N).
179*> \endverbatim
180*>
181*> \param[in,out] EQUED
182*> \verbatim
183*> EQUED is CHARACTER*1
184*> Specifies the form of equilibration that was done.
185*> = 'N': No equilibration (always true if FACT = 'N').
186*> = 'Y': Equilibration was done, i.e., A has been replaced by
187*> diag(S) * A * diag(S).
188*> EQUED is an input argument if FACT = 'F'; otherwise, it is an
189*> output argument.
190*> \endverbatim
191*>
192*> \param[in,out] S
193*> \verbatim
194*> S is DOUBLE PRECISION array, dimension (N)
195*> The scale factors for A; not accessed if EQUED = 'N'. S is
196*> an input argument if FACT = 'F'; otherwise, S is an output
197*> argument. If FACT = 'F' and EQUED = 'Y', each element of S
198*> must be positive.
199*> \endverbatim
200*>
201*> \param[in,out] B
202*> \verbatim
203*> B is COMPLEX*16 array, dimension (LDB,NRHS)
204*> On entry, the N-by-NRHS righthand side matrix B.
205*> On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y',
206*> B is overwritten by diag(S) * B.
207*> \endverbatim
208*>
209*> \param[in] LDB
210*> \verbatim
211*> LDB is INTEGER
212*> The leading dimension of the array B. LDB >= max(1,N).
213*> \endverbatim
214*>
215*> \param[out] X
216*> \verbatim
217*> X is COMPLEX*16 array, dimension (LDX,NRHS)
218*> If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X to
219*> the original system of equations. Note that if EQUED = 'Y',
220*> A and B are modified on exit, and the solution to the
221*> equilibrated system is inv(diag(S))*X.
222*> \endverbatim
223*>
224*> \param[in] LDX
225*> \verbatim
226*> LDX is INTEGER
227*> The leading dimension of the array X. LDX >= max(1,N).
228*> \endverbatim
229*>
230*> \param[out] RCOND
231*> \verbatim
232*> RCOND is DOUBLE PRECISION
233*> The estimate of the reciprocal condition number of the matrix
234*> A after equilibration (if done). If RCOND is less than the
235*> machine precision (in particular, if RCOND = 0), the matrix
236*> is singular to working precision. This condition is
237*> indicated by a return code of INFO > 0.
238*> \endverbatim
239*>
240*> \param[out] FERR
241*> \verbatim
242*> FERR is DOUBLE PRECISION array, dimension (NRHS)
243*> The estimated forward error bound for each solution vector
244*> X(j) (the j-th column of the solution matrix X).
245*> If XTRUE is the true solution corresponding to X(j), FERR(j)
246*> is an estimated upper bound for the magnitude of the largest
247*> element in (X(j) - XTRUE) divided by the magnitude of the
248*> largest element in X(j). The estimate is as reliable as
249*> the estimate for RCOND, and is almost always a slight
250*> overestimate of the true error.
251*> \endverbatim
252*>
253*> \param[out] BERR
254*> \verbatim
255*> BERR is DOUBLE PRECISION array, dimension (NRHS)
256*> The componentwise relative backward error of each solution
257*> vector X(j) (i.e., the smallest relative change in
258*> any element of A or B that makes X(j) an exact solution).
259*> \endverbatim
260*>
261*> \param[out] WORK
262*> \verbatim
263*> WORK is COMPLEX*16 array, dimension (2*N)
264*> \endverbatim
265*>
266*> \param[out] RWORK
267*> \verbatim
268*> RWORK is DOUBLE PRECISION array, dimension (N)
269*> \endverbatim
270*>
271*> \param[out] INFO
272*> \verbatim
273*> INFO is INTEGER
274*> = 0: successful exit
275*> < 0: if INFO = -i, the i-th argument had an illegal value
276*> > 0: if INFO = i, and i is
277*> <= N: the leading principal minor of order i of A
278*> is not positive, so the factorization could not
279*> be completed, and the solution has not been
280*> computed. RCOND = 0 is returned.
281*> = N+1: U is nonsingular, but RCOND is less than machine
282*> precision, meaning that the matrix is singular
283*> to working precision. Nevertheless, the
284*> solution and error bounds are computed because
285*> there are a number of situations where the
286*> computed solution can be more accurate than the
287*> value of RCOND would suggest.
288*> \endverbatim
289*
290* Authors:
291* ========
292*
293*> \author Univ. of Tennessee
294*> \author Univ. of California Berkeley
295*> \author Univ. of Colorado Denver
296*> \author NAG Ltd.
297*
298*> \ingroup posvx
299*
300* =====================================================================
301 SUBROUTINE zposvx( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF,
302 $ EQUED,
303 $ S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK,
304 $ RWORK, INFO )
305*
306* -- LAPACK driver routine --
307* -- LAPACK is a software package provided by Univ. of Tennessee, --
308* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
309*
310* .. Scalar Arguments ..
311 CHARACTER EQUED, FACT, UPLO
312 INTEGER INFO, LDA, LDAF, LDB, LDX, N, NRHS
313 DOUBLE PRECISION RCOND
314* ..
315* .. Array Arguments ..
316 DOUBLE PRECISION BERR( * ), FERR( * ), RWORK( * ), S( * )
317 COMPLEX*16 A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
318 $ WORK( * ), X( LDX, * )
319* ..
320*
321* =====================================================================
322*
323* .. Parameters ..
324 DOUBLE PRECISION ZERO, ONE
325 PARAMETER ( ZERO = 0.0d+0, one = 1.0d+0 )
326* ..
327* .. Local Scalars ..
328 LOGICAL EQUIL, NOFACT, RCEQU
329 INTEGER I, INFEQU, J
330 DOUBLE PRECISION AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM
331* ..
332* .. External Functions ..
333 LOGICAL LSAME
334 DOUBLE PRECISION DLAMCH, ZLANHE
335 EXTERNAL LSAME, DLAMCH, ZLANHE
336* ..
337* .. External Subroutines ..
338 EXTERNAL xerbla, zlacpy, zlaqhe, zpocon, zpoequ,
339 $ zporfs,
340 $ zpotrf, zpotrs
341* ..
342* .. Intrinsic Functions ..
343 INTRINSIC max, min
344* ..
345* .. Executable Statements ..
346*
347 info = 0
348 nofact = lsame( fact, 'N' )
349 equil = lsame( fact, 'E' )
350 IF( nofact .OR. equil ) THEN
351 equed = 'N'
352 rcequ = .false.
353 ELSE
354 rcequ = lsame( equed, 'Y' )
355 smlnum = dlamch( 'Safe minimum' )
356 bignum = one / smlnum
357 END IF
358*
359* Test the input parameters.
360*
361 IF( .NOT.nofact .AND.
362 $ .NOT.equil .AND.
363 $ .NOT.lsame( fact, 'F' ) )
364 $ THEN
365 info = -1
366 ELSE IF( .NOT.lsame( uplo, 'U' ) .AND.
367 $ .NOT.lsame( uplo, 'L' ) )
368 $ THEN
369 info = -2
370 ELSE IF( n.LT.0 ) THEN
371 info = -3
372 ELSE IF( nrhs.LT.0 ) THEN
373 info = -4
374 ELSE IF( lda.LT.max( 1, n ) ) THEN
375 info = -6
376 ELSE IF( ldaf.LT.max( 1, n ) ) THEN
377 info = -8
378 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
379 $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
380 info = -9
381 ELSE
382 IF( rcequ ) THEN
383 smin = bignum
384 smax = zero
385 DO 10 j = 1, n
386 smin = min( smin, s( j ) )
387 smax = max( smax, s( j ) )
388 10 CONTINUE
389 IF( smin.LE.zero ) THEN
390 info = -10
391 ELSE IF( n.GT.0 ) THEN
392 scond = max( smin, smlnum ) / min( smax, bignum )
393 ELSE
394 scond = one
395 END IF
396 END IF
397 IF( info.EQ.0 ) THEN
398 IF( ldb.LT.max( 1, n ) ) THEN
399 info = -12
400 ELSE IF( ldx.LT.max( 1, n ) ) THEN
401 info = -14
402 END IF
403 END IF
404 END IF
405*
406 IF( info.NE.0 ) THEN
407 CALL xerbla( 'ZPOSVX', -info )
408 RETURN
409 END IF
410*
411 IF( equil ) THEN
412*
413* Compute row and column scalings to equilibrate the matrix A.
414*
415 CALL zpoequ( n, a, lda, s, scond, amax, infequ )
416 IF( infequ.EQ.0 ) THEN
417*
418* Equilibrate the matrix.
419*
420 CALL zlaqhe( uplo, n, a, lda, s, scond, amax, equed )
421 rcequ = lsame( equed, 'Y' )
422 END IF
423 END IF
424*
425* Scale the right hand side.
426*
427 IF( rcequ ) THEN
428 DO 30 j = 1, nrhs
429 DO 20 i = 1, n
430 b( i, j ) = s( i )*b( i, j )
431 20 CONTINUE
432 30 CONTINUE
433 END IF
434*
435 IF( nofact .OR. equil ) THEN
436*
437* Compute the Cholesky factorization A = U**H *U or A = L*L**H.
438*
439 CALL zlacpy( uplo, n, n, a, lda, af, ldaf )
440 CALL zpotrf( uplo, n, af, ldaf, info )
441*
442* Return if INFO is non-zero.
443*
444 IF( info.GT.0 )THEN
445 rcond = zero
446 RETURN
447 END IF
448 END IF
449*
450* Compute the norm of the matrix A.
451*
452 anorm = zlanhe( '1', uplo, n, a, lda, rwork )
453*
454* Compute the reciprocal of the condition number of A.
455*
456 CALL zpocon( uplo, n, af, ldaf, anorm, rcond, work, rwork,
457 $ info )
458*
459* Compute the solution matrix X.
460*
461 CALL zlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
462 CALL zpotrs( uplo, n, nrhs, af, ldaf, x, ldx, info )
463*
464* Use iterative refinement to improve the computed solution and
465* compute error bounds and backward error estimates for it.
466*
467 CALL zporfs( uplo, n, nrhs, a, lda, af, ldaf, b, ldb, x, ldx,
468 $ ferr, berr, work, rwork, info )
469*
470* Transform the solution matrix X to a solution of the original
471* system.
472*
473 IF( rcequ ) THEN
474 DO 50 j = 1, nrhs
475 DO 40 i = 1, n
476 x( i, j ) = s( i )*x( i, j )
477 40 CONTINUE
478 50 CONTINUE
479 DO 60 j = 1, nrhs
480 ferr( j ) = ferr( j ) / scond
481 60 CONTINUE
482 END IF
483*
484* Set INFO = N+1 if the matrix is singular to working precision.
485*
486 IF( rcond.LT.dlamch( 'Epsilon' ) )
487 $ info = n + 1
488*
489 RETURN
490*
491* End of ZPOSVX
492*
493 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:101
subroutine zlaqhe(uplo, n, a, lda, s, scond, amax, equed)
ZLAQHE scales a Hermitian matrix.
Definition zlaqhe.f:132
subroutine zpocon(uplo, n, a, lda, anorm, rcond, work, rwork, info)
ZPOCON
Definition zpocon.f:119
subroutine zpoequ(n, a, lda, s, scond, amax, info)
ZPOEQU
Definition zpoequ.f:111
subroutine zporfs(uplo, n, nrhs, a, lda, af, ldaf, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZPORFS
Definition zporfs.f:181
subroutine zposvx(fact, uplo, n, nrhs, a, lda, af, ldaf, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
ZPOSVX computes the solution to system of linear equations A * X = B for PO matrices
Definition zposvx.f:305
subroutine zpotrf(uplo, n, a, lda, info)
ZPOTRF
Definition zpotrf.f:105
subroutine zpotrs(uplo, n, nrhs, a, lda, b, ldb, info)
ZPOTRS
Definition zpotrs.f:108