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