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