b653584b91627bcb8129afd4d89f0e586153e7f2
[platform/upstream/lapack.git] / SRC / cppsvx.f
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 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 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 minor of order i of A is
266 *>                       not positive definite, so the factorization
267 *>                       could not be completed, and the solution has not
268 *>                       been 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 *> \date April 2012
287 *
288 *> \ingroup complexOTHERsolve
289 *
290 *> \par Further Details:
291 *  =====================
292 *>
293 *> \verbatim
294 *>
295 *>  The packed storage scheme is illustrated by the following example
296 *>  when N = 4, UPLO = 'U':
297 *>
298 *>  Two-dimensional storage of the Hermitian matrix A:
299 *>
300 *>     a11 a12 a13 a14
301 *>         a22 a23 a24
302 *>             a33 a34     (aij = conjg(aji))
303 *>                 a44
304 *>
305 *>  Packed storage of the upper triangle of A:
306 *>
307 *>  AP = [ a11, a12, a22, a13, a23, a33, a14, a24, a34, a44 ]
308 *> \endverbatim
309 *>
310 *  =====================================================================
311       SUBROUTINE CPPSVX( FACT, UPLO, N, NRHS, AP, AFP, EQUED, S, B, LDB,
312      $                   X, LDX, RCOND, FERR, BERR, WORK, RWORK, INFO )
313 *
314 *  -- LAPACK driver routine (version 3.4.1) --
315 *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
316 *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
317 *     April 2012
318 *
319 *     .. Scalar Arguments ..
320       CHARACTER          EQUED, FACT, UPLO
321       INTEGER            INFO, LDB, LDX, N, NRHS
322       REAL               RCOND
323 *     ..
324 *     .. Array Arguments ..
325       REAL               BERR( * ), FERR( * ), RWORK( * ), S( * )
326       COMPLEX            AFP( * ), AP( * ), B( LDB, * ), WORK( * ),
327      $                   X( LDX, * )
328 *     ..
329 *
330 *  =====================================================================
331 *
332 *     .. Parameters ..
333       REAL               ZERO, ONE
334       PARAMETER          ( ZERO = 0.0E+0, ONE = 1.0E+0 )
335 *     ..
336 *     .. Local Scalars ..
337       LOGICAL            EQUIL, NOFACT, RCEQU
338       INTEGER            I, INFEQU, J
339       REAL               AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM
340 *     ..
341 *     .. External Functions ..
342       LOGICAL            LSAME
343       REAL               CLANHP, SLAMCH
344       EXTERNAL           LSAME, CLANHP, SLAMCH
345 *     ..
346 *     .. External Subroutines ..
347       EXTERNAL           CCOPY, CLACPY, CLAQHP, CPPCON, CPPEQU, CPPRFS,
348      $                   CPPTRF, CPPTRS, XERBLA
349 *     ..
350 *     .. Intrinsic Functions ..
351       INTRINSIC          MAX, MIN
352 *     ..
353 *     .. Executable Statements ..
354 *
355       INFO = 0
356       NOFACT = LSAME( FACT, 'N' )
357       EQUIL = LSAME( FACT, 'E' )
358       IF( NOFACT .OR. EQUIL ) THEN
359          EQUED = 'N'
360          RCEQU = .FALSE.
361       ELSE
362          RCEQU = LSAME( EQUED, 'Y' )
363          SMLNUM = SLAMCH( 'Safe minimum' )
364          BIGNUM = ONE / SMLNUM
365       END IF
366 *
367 *     Test the input parameters.
368 *
369       IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.LSAME( FACT, 'F' ) )
370      $     THEN
371          INFO = -1
372       ELSE IF( .NOT.LSAME( UPLO, 'U' ) .AND. .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( 'CPPSVX', -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 CPPEQU( UPLO, N, AP, S, SCOND, AMAX, INFEQU )
417          IF( INFEQU.EQ.0 ) THEN
418 *
419 *           Equilibrate the matrix.
420 *
421             CALL CLAQHP( 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**H * U or A = L * L**H.
439 *
440          CALL CCOPY( N*( N+1 ) / 2, AP, 1, AFP, 1 )
441          CALL CPPTRF( 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 = CLANHP( 'I', UPLO, N, AP, RWORK )
454 *
455 *     Compute the reciprocal of the condition number of A.
456 *
457       CALL CPPCON( UPLO, N, AFP, ANORM, RCOND, WORK, RWORK, INFO )
458 *
459 *     Compute the solution matrix X.
460 *
461       CALL CLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
462       CALL CPPTRS( 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 CPPRFS( UPLO, N, NRHS, AP, AFP, B, LDB, X, LDX, FERR, BERR,
468      $             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.SLAMCH( 'Epsilon' ) )
487      $   INFO = N + 1
488 *
489       RETURN
490 *
491 *     End of CPPSVX
492 *
493       END