24144f38ab8197073d244987da8bd3420479c403
[platform/upstream/lapack.git] / SRC / sppsvx.f
1 *> \brief <b> SPPSVX 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 SPPSVX + dependencies 
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sppsvx.f"> 
11 *> [TGZ]</a> 
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sppsvx.f"> 
13 *> [ZIP]</a> 
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sppsvx.f"> 
15 *> [TXT]</a>
16 *> \endhtmlonly 
17 *
18 *  Definition:
19 *  ===========
20 *
21 *       SUBROUTINE SPPSVX( FACT, UPLO, N, NRHS, AP, AFP, EQUED, S, B, LDB,
22 *                          X, LDX, RCOND, FERR, BERR, WORK, IWORK, 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 *       INTEGER            IWORK( * )
31 *       REAL               AFP( * ), AP( * ), B( LDB, * ), BERR( * ),
32 *      $                   FERR( * ), S( * ), WORK( * ), X( LDX, * )
33 *       ..
34 *  
35 *
36 *> \par Purpose:
37 *  =============
38 *>
39 *> \verbatim
40 *>
41 *> SPPSVX uses the Cholesky factorization A = U**T*U or A = L*L**T to
42 *> compute the solution to a real system of linear equations
43 *>    A * X = B,
44 *> where A is an N-by-N symmetric 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**T* U,  if UPLO = 'U', or
68 *>       A = L * L**T,  if UPLO = 'L',
69 *>    where U is an upper triangular matrix and L is a lower triangular
70 *>    matrix.
71 *>
72 *> 3. If the leading i-by-i principal minor is not positive definite,
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 REAL array, dimension (N*(N+1)/2)
134 *>          On entry, the upper or lower triangle of the symmetric 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 REAL array, dimension
151 *>                            (N*(N+1)/2)
152 *>          If FACT = 'F', then AFP is an input argument and on entry
153 *>          contains the triangular factor U or L from the Cholesky
154 *>          factorization A = U**T*U or A = L*L**T, in the same storage
155 *>          format as A.  If EQUED .ne. 'N', then AFP is the factored
156 *>          form of the equilibrated matrix A.
157 *>
158 *>          If FACT = 'N', then AFP is an output argument and on exit
159 *>          returns the triangular factor U or L from the Cholesky
160 *>          factorization A = U**T * U or A = L * L**T of the original
161 *>          matrix A.
162 *>
163 *>          If FACT = 'E', then AFP is an output argument and on exit
164 *>          returns the triangular factor U or L from the Cholesky
165 *>          factorization A = U**T * U or A = L * L**T of the equilibrated
166 *>          matrix A (see the description of AP for the form of the
167 *>          equilibrated matrix).
168 *> \endverbatim
169 *>
170 *> \param[in,out] EQUED
171 *> \verbatim
172 *>          EQUED is CHARACTER*1
173 *>          Specifies the form of equilibration that was done.
174 *>          = 'N':  No equilibration (always true if FACT = 'N').
175 *>          = 'Y':  Equilibration was done, i.e., A has been replaced by
176 *>                  diag(S) * A * diag(S).
177 *>          EQUED is an input argument if FACT = 'F'; otherwise, it is an
178 *>          output argument.
179 *> \endverbatim
180 *>
181 *> \param[in,out] S
182 *> \verbatim
183 *>          S is REAL array, dimension (N)
184 *>          The scale factors for A; not accessed if EQUED = 'N'.  S is
185 *>          an input argument if FACT = 'F'; otherwise, S is an output
186 *>          argument.  If FACT = 'F' and EQUED = 'Y', each element of S
187 *>          must be positive.
188 *> \endverbatim
189 *>
190 *> \param[in,out] B
191 *> \verbatim
192 *>          B is REAL array, dimension (LDB,NRHS)
193 *>          On entry, the N-by-NRHS right hand side matrix B.
194 *>          On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y',
195 *>          B is overwritten by diag(S) * B.
196 *> \endverbatim
197 *>
198 *> \param[in] LDB
199 *> \verbatim
200 *>          LDB is INTEGER
201 *>          The leading dimension of the array B.  LDB >= max(1,N).
202 *> \endverbatim
203 *>
204 *> \param[out] X
205 *> \verbatim
206 *>          X is REAL array, dimension (LDX,NRHS)
207 *>          If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X to
208 *>          the original system of equations.  Note that if EQUED = 'Y',
209 *>          A and B are modified on exit, and the solution to the
210 *>          equilibrated system is inv(diag(S))*X.
211 *> \endverbatim
212 *>
213 *> \param[in] LDX
214 *> \verbatim
215 *>          LDX is INTEGER
216 *>          The leading dimension of the array X.  LDX >= max(1,N).
217 *> \endverbatim
218 *>
219 *> \param[out] RCOND
220 *> \verbatim
221 *>          RCOND is REAL
222 *>          The estimate of the reciprocal condition number of the matrix
223 *>          A after equilibration (if done).  If RCOND is less than the
224 *>          machine precision (in particular, if RCOND = 0), the matrix
225 *>          is singular to working precision.  This condition is
226 *>          indicated by a return code of INFO > 0.
227 *> \endverbatim
228 *>
229 *> \param[out] FERR
230 *> \verbatim
231 *>          FERR is REAL array, dimension (NRHS)
232 *>          The estimated forward error bound for each solution vector
233 *>          X(j) (the j-th column of the solution matrix X).
234 *>          If XTRUE is the true solution corresponding to X(j), FERR(j)
235 *>          is an estimated upper bound for the magnitude of the largest
236 *>          element in (X(j) - XTRUE) divided by the magnitude of the
237 *>          largest element in X(j).  The estimate is as reliable as
238 *>          the estimate for RCOND, and is almost always a slight
239 *>          overestimate of the true error.
240 *> \endverbatim
241 *>
242 *> \param[out] BERR
243 *> \verbatim
244 *>          BERR is REAL array, dimension (NRHS)
245 *>          The componentwise relative backward error of each solution
246 *>          vector X(j) (i.e., the smallest relative change in
247 *>          any element of A or B that makes X(j) an exact solution).
248 *> \endverbatim
249 *>
250 *> \param[out] WORK
251 *> \verbatim
252 *>          WORK is REAL array, dimension (3*N)
253 *> \endverbatim
254 *>
255 *> \param[out] IWORK
256 *> \verbatim
257 *>          IWORK is INTEGER array, dimension (N)
258 *> \endverbatim
259 *>
260 *> \param[out] INFO
261 *> \verbatim
262 *>          INFO is INTEGER
263 *>          = 0:  successful exit
264 *>          < 0:  if INFO = -i, the i-th argument had an illegal value
265 *>          > 0:  if INFO = i, and i is
266 *>                <= N:  the leading minor of order i of A is
267 *>                       not positive definite, so the factorization
268 *>                       could not be completed, and the solution has not
269 *>                       been computed. RCOND = 0 is returned.
270 *>                = N+1: U is nonsingular, but RCOND is less than machine
271 *>                       precision, meaning that the matrix is singular
272 *>                       to working precision.  Nevertheless, the
273 *>                       solution and error bounds are computed because
274 *>                       there are a number of situations where the
275 *>                       computed solution can be more accurate than the
276 *>                       value of RCOND would suggest.
277 *> \endverbatim
278 *
279 *  Authors:
280 *  ========
281 *
282 *> \author Univ. of Tennessee 
283 *> \author Univ. of California Berkeley 
284 *> \author Univ. of Colorado Denver 
285 *> \author NAG Ltd. 
286 *
287 *> \date April 2012
288 *
289 *> \ingroup realOTHERsolve
290 *
291 *> \par Further Details:
292 *  =====================
293 *>
294 *> \verbatim
295 *>
296 *>  The packed storage scheme is illustrated by the following example
297 *>  when N = 4, UPLO = 'U':
298 *>
299 *>  Two-dimensional storage of the symmetric matrix A:
300 *>
301 *>     a11 a12 a13 a14
302 *>         a22 a23 a24
303 *>             a33 a34     (aij = conjg(aji))
304 *>                 a44
305 *>
306 *>  Packed storage of the upper triangle of A:
307 *>
308 *>  AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]
309 *> \endverbatim
310 *>
311 *  =====================================================================
312       SUBROUTINE SPPSVX( FACT, UPLO, N, NRHS, AP, AFP, EQUED, S, B, LDB,
313      $                   X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO )
314 *
315 *  -- LAPACK driver routine (version 3.4.1) --
316 *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
317 *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
318 *     April 2012
319 *
320 *     .. Scalar Arguments ..
321       CHARACTER          EQUED, FACT, UPLO
322       INTEGER            INFO, LDB, LDX, N, NRHS
323       REAL               RCOND
324 *     ..
325 *     .. Array Arguments ..
326       INTEGER            IWORK( * )
327       REAL               AFP( * ), AP( * ), B( LDB, * ), BERR( * ),
328      $                   FERR( * ), S( * ), WORK( * ), X( LDX, * )
329 *     ..
330 *
331 *  =====================================================================
332 *
333 *     .. Parameters ..
334       REAL               ZERO, ONE
335       PARAMETER          ( ZERO = 0.0E+0, ONE = 1.0E+0 )
336 *     ..
337 *     .. Local Scalars ..
338       LOGICAL            EQUIL, NOFACT, RCEQU
339       INTEGER            I, INFEQU, J
340       REAL               AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM
341 *     ..
342 *     .. External Functions ..
343       LOGICAL            LSAME
344       REAL               SLAMCH, SLANSP
345       EXTERNAL           LSAME, SLAMCH, SLANSP
346 *     ..
347 *     .. External Subroutines ..
348       EXTERNAL           SCOPY, SLACPY, SLAQSP, SPPCON, SPPEQU, SPPRFS,
349      $                   SPPTRF, SPPTRS, XERBLA
350 *     ..
351 *     .. Intrinsic Functions ..
352       INTRINSIC          MAX, MIN
353 *     ..
354 *     .. Executable Statements ..
355 *
356       INFO = 0
357       NOFACT = LSAME( FACT, 'N' )
358       EQUIL = LSAME( FACT, 'E' )
359       IF( NOFACT .OR. EQUIL ) THEN
360          EQUED = 'N'
361          RCEQU = .FALSE.
362       ELSE
363          RCEQU = LSAME( EQUED, 'Y' )
364          SMLNUM = SLAMCH( 'Safe minimum' )
365          BIGNUM = ONE / SMLNUM
366       END IF
367 *
368 *     Test the input parameters.
369 *
370       IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.LSAME( FACT, 'F' ) )
371      $     THEN
372          INFO = -1
373       ELSE IF( .NOT.LSAME( UPLO, 'U' ) .AND. .NOT.LSAME( UPLO, 'L' ) )
374      $          THEN
375          INFO = -2
376       ELSE IF( N.LT.0 ) THEN
377          INFO = -3
378       ELSE IF( NRHS.LT.0 ) THEN
379          INFO = -4
380       ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
381      $         ( RCEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
382          INFO = -7
383       ELSE
384          IF( RCEQU ) THEN
385             SMIN = BIGNUM
386             SMAX = ZERO
387             DO 10 J = 1, N
388                SMIN = MIN( SMIN, S( J ) )
389                SMAX = MAX( SMAX, S( J ) )
390    10       CONTINUE
391             IF( SMIN.LE.ZERO ) THEN
392                INFO = -8
393             ELSE IF( N.GT.0 ) THEN
394                SCOND = MAX( SMIN, SMLNUM ) / MIN( SMAX, BIGNUM )
395             ELSE
396                SCOND = ONE
397             END IF
398          END IF
399          IF( INFO.EQ.0 ) THEN
400             IF( LDB.LT.MAX( 1, N ) ) THEN
401                INFO = -10
402             ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
403                INFO = -12
404             END IF
405          END IF
406       END IF
407 *
408       IF( INFO.NE.0 ) THEN
409          CALL XERBLA( 'SPPSVX', -INFO )
410          RETURN
411       END IF
412 *
413       IF( EQUIL ) THEN
414 *
415 *        Compute row and column scalings to equilibrate the matrix A.
416 *
417          CALL SPPEQU( UPLO, N, AP, S, SCOND, AMAX, INFEQU )
418          IF( INFEQU.EQ.0 ) THEN
419 *
420 *           Equilibrate the matrix.
421 *
422             CALL SLAQSP( UPLO, N, AP, S, SCOND, AMAX, EQUED )
423             RCEQU = LSAME( EQUED, 'Y' )
424          END IF
425       END IF
426 *
427 *     Scale the right-hand side.
428 *
429       IF( RCEQU ) THEN
430          DO 30 J = 1, NRHS
431             DO 20 I = 1, N
432                B( I, J ) = S( I )*B( I, J )
433    20       CONTINUE
434    30    CONTINUE
435       END IF
436 *
437       IF( NOFACT .OR. EQUIL ) THEN
438 *
439 *        Compute the Cholesky factorization A = U**T * U or A = L * L**T.
440 *
441          CALL SCOPY( N*( N+1 ) / 2, AP, 1, AFP, 1 )
442          CALL SPPTRF( UPLO, N, AFP, INFO )
443 *
444 *        Return if INFO is non-zero.
445 *
446          IF( INFO.GT.0 )THEN
447             RCOND = ZERO
448             RETURN
449          END IF
450       END IF
451 *
452 *     Compute the norm of the matrix A.
453 *
454       ANORM = SLANSP( 'I', UPLO, N, AP, WORK )
455 *
456 *     Compute the reciprocal of the condition number of A.
457 *
458       CALL SPPCON( UPLO, N, AFP, ANORM, RCOND, WORK, IWORK, INFO )
459 *
460 *     Compute the solution matrix X.
461 *
462       CALL SLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
463       CALL SPPTRS( UPLO, N, NRHS, AFP, X, LDX, INFO )
464 *
465 *     Use iterative refinement to improve the computed solution and
466 *     compute error bounds and backward error estimates for it.
467 *
468       CALL SPPRFS( UPLO, N, NRHS, AP, AFP, B, LDB, X, LDX, FERR, 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.SLAMCH( 'Epsilon' ) )
488      $   INFO = N + 1
489 *
490       RETURN
491 *
492 *     End of SPPSVX
493 *
494       END