C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / ssprfs.c
1 #include <math.h>
2 #include <stdlib.h>
3 #include <string.h>
4 #include <stdio.h>
5 #include <complex.h>
6 #ifdef complex
7 #undef complex
8 #endif
9 #ifdef I
10 #undef I
11 #endif
12
13 #if defined(_WIN64)
14 typedef long long BLASLONG;
15 typedef unsigned long long BLASULONG;
16 #else
17 typedef long BLASLONG;
18 typedef unsigned long BLASULONG;
19 #endif
20
21 #ifdef LAPACK_ILP64
22 typedef BLASLONG blasint;
23 #if defined(_WIN64)
24 #define blasabs(x) llabs(x)
25 #else
26 #define blasabs(x) labs(x)
27 #endif
28 #else
29 typedef int blasint;
30 #define blasabs(x) abs(x)
31 #endif
32
33 typedef blasint integer;
34
35 typedef unsigned int uinteger;
36 typedef char *address;
37 typedef short int shortint;
38 typedef float real;
39 typedef double doublereal;
40 typedef struct { real r, i; } complex;
41 typedef struct { doublereal r, i; } doublecomplex;
42 #ifdef _MSC_VER
43 static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
44 static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
45 static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
46 static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
47 #else
48 static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
49 static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
50 static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
51 static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
52 #endif
53 #define pCf(z) (*_pCf(z))
54 #define pCd(z) (*_pCd(z))
55 typedef int logical;
56 typedef short int shortlogical;
57 typedef char logical1;
58 typedef char integer1;
59
60 #define TRUE_ (1)
61 #define FALSE_ (0)
62
63 /* Extern is for use with -E */
64 #ifndef Extern
65 #define Extern extern
66 #endif
67
68 /* I/O stuff */
69
70 typedef int flag;
71 typedef int ftnlen;
72 typedef int ftnint;
73
74 /*external read, write*/
75 typedef struct
76 {       flag cierr;
77         ftnint ciunit;
78         flag ciend;
79         char *cifmt;
80         ftnint cirec;
81 } cilist;
82
83 /*internal read, write*/
84 typedef struct
85 {       flag icierr;
86         char *iciunit;
87         flag iciend;
88         char *icifmt;
89         ftnint icirlen;
90         ftnint icirnum;
91 } icilist;
92
93 /*open*/
94 typedef struct
95 {       flag oerr;
96         ftnint ounit;
97         char *ofnm;
98         ftnlen ofnmlen;
99         char *osta;
100         char *oacc;
101         char *ofm;
102         ftnint orl;
103         char *oblnk;
104 } olist;
105
106 /*close*/
107 typedef struct
108 {       flag cerr;
109         ftnint cunit;
110         char *csta;
111 } cllist;
112
113 /*rewind, backspace, endfile*/
114 typedef struct
115 {       flag aerr;
116         ftnint aunit;
117 } alist;
118
119 /* inquire */
120 typedef struct
121 {       flag inerr;
122         ftnint inunit;
123         char *infile;
124         ftnlen infilen;
125         ftnint  *inex;  /*parameters in standard's order*/
126         ftnint  *inopen;
127         ftnint  *innum;
128         ftnint  *innamed;
129         char    *inname;
130         ftnlen  innamlen;
131         char    *inacc;
132         ftnlen  inacclen;
133         char    *inseq;
134         ftnlen  inseqlen;
135         char    *indir;
136         ftnlen  indirlen;
137         char    *infmt;
138         ftnlen  infmtlen;
139         char    *inform;
140         ftnint  informlen;
141         char    *inunf;
142         ftnlen  inunflen;
143         ftnint  *inrecl;
144         ftnint  *innrec;
145         char    *inblank;
146         ftnlen  inblanklen;
147 } inlist;
148
149 #define VOID void
150
151 union Multitype {       /* for multiple entry points */
152         integer1 g;
153         shortint h;
154         integer i;
155         /* longint j; */
156         real r;
157         doublereal d;
158         complex c;
159         doublecomplex z;
160         };
161
162 typedef union Multitype Multitype;
163
164 struct Vardesc {        /* for Namelist */
165         char *name;
166         char *addr;
167         ftnlen *dims;
168         int  type;
169         };
170 typedef struct Vardesc Vardesc;
171
172 struct Namelist {
173         char *name;
174         Vardesc **vars;
175         int nvars;
176         };
177 typedef struct Namelist Namelist;
178
179 #define abs(x) ((x) >= 0 ? (x) : -(x))
180 #define dabs(x) (fabs(x))
181 #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
182 #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
183 #define dmin(a,b) (f2cmin(a,b))
184 #define dmax(a,b) (f2cmax(a,b))
185 #define bit_test(a,b)   ((a) >> (b) & 1)
186 #define bit_clear(a,b)  ((a) & ~((uinteger)1 << (b)))
187 #define bit_set(a,b)    ((a) |  ((uinteger)1 << (b)))
188
189 #define abort_() { sig_die("Fortran abort routine called", 1); }
190 #define c_abs(z) (cabsf(Cf(z)))
191 #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
192 #ifdef _MSC_VER
193 #define c_div(c, a, b) {Cf(c)._Val[0] = (Cf(a)._Val[0]/Cf(b)._Val[0]); Cf(c)._Val[1]=(Cf(a)._Val[1]/Cf(b)._Val[1]);}
194 #define z_div(c, a, b) {Cd(c)._Val[0] = (Cd(a)._Val[0]/Cd(b)._Val[0]); Cd(c)._Val[1]=(Cd(a)._Val[1]/df(b)._Val[1]);}
195 #else
196 #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
197 #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
198 #endif
199 #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
200 #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
201 #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
202 //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
203 #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
204 #define d_abs(x) (fabs(*(x)))
205 #define d_acos(x) (acos(*(x)))
206 #define d_asin(x) (asin(*(x)))
207 #define d_atan(x) (atan(*(x)))
208 #define d_atn2(x, y) (atan2(*(x),*(y)))
209 #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
210 #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
211 #define d_cos(x) (cos(*(x)))
212 #define d_cosh(x) (cosh(*(x)))
213 #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
214 #define d_exp(x) (exp(*(x)))
215 #define d_imag(z) (cimag(Cd(z)))
216 #define r_imag(z) (cimagf(Cf(z)))
217 #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
218 #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
219 #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
220 #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
221 #define d_log(x) (log(*(x)))
222 #define d_mod(x, y) (fmod(*(x), *(y)))
223 #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
224 #define d_nint(x) u_nint(*(x))
225 #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
226 #define d_sign(a,b) u_sign(*(a),*(b))
227 #define r_sign(a,b) u_sign(*(a),*(b))
228 #define d_sin(x) (sin(*(x)))
229 #define d_sinh(x) (sinh(*(x)))
230 #define d_sqrt(x) (sqrt(*(x)))
231 #define d_tan(x) (tan(*(x)))
232 #define d_tanh(x) (tanh(*(x)))
233 #define i_abs(x) abs(*(x))
234 #define i_dnnt(x) ((integer)u_nint(*(x)))
235 #define i_len(s, n) (n)
236 #define i_nint(x) ((integer)u_nint(*(x)))
237 #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
238 #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
239 #define pow_si(B,E) spow_ui(*(B),*(E))
240 #define pow_ri(B,E) spow_ui(*(B),*(E))
241 #define pow_di(B,E) dpow_ui(*(B),*(E))
242 #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
243 #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
244 #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
245 #define s_cat(lpp, rpp, rnp, np, llp) {         ftnlen i, nc, ll; char *f__rp, *lp;     ll = (llp); lp = (lpp);         for(i=0; i < (int)*(np); ++i) {                 nc = ll;                if((rnp)[i] < nc) nc = (rnp)[i];                ll -= nc;               f__rp = (rpp)[i];               while(--nc >= 0) *lp++ = *(f__rp)++;         }  while(--ll >= 0) *lp++ = ' '; }
246 #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
247 #define s_copy(A,B,C,D) { int __i,__m; for (__i=0, __m=f2cmin((C),(D)); __i<__m && (B)[__i] != 0; ++__i) (A)[__i] = (B)[__i]; }
248 #define sig_die(s, kill) { exit(1); }
249 #define s_stop(s, n) {exit(0);}
250 static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
251 #define z_abs(z) (cabs(Cd(z)))
252 #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
253 #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
254 #define myexit_() break;
255 #define mycycle() continue;
256 #define myceiling(w) {ceil(w)}
257 #define myhuge(w) {HUGE_VAL}
258 //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
259 #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
260
261 /* procedure parameter types for -A and -C++ */
262
263 #define F2C_proc_par_types 1
264 #ifdef __cplusplus
265 typedef logical (*L_fp)(...);
266 #else
267 typedef logical (*L_fp)();
268 #endif
269
270 static float spow_ui(float x, integer n) {
271         float pow=1.0; unsigned long int u;
272         if(n != 0) {
273                 if(n < 0) n = -n, x = 1/x;
274                 for(u = n; ; ) {
275                         if(u & 01) pow *= x;
276                         if(u >>= 1) x *= x;
277                         else break;
278                 }
279         }
280         return pow;
281 }
282 static double dpow_ui(double x, integer n) {
283         double pow=1.0; unsigned long int u;
284         if(n != 0) {
285                 if(n < 0) n = -n, x = 1/x;
286                 for(u = n; ; ) {
287                         if(u & 01) pow *= x;
288                         if(u >>= 1) x *= x;
289                         else break;
290                 }
291         }
292         return pow;
293 }
294 #ifdef _MSC_VER
295 static _Fcomplex cpow_ui(complex x, integer n) {
296         complex pow={1.0,0.0}; unsigned long int u;
297                 if(n != 0) {
298                 if(n < 0) n = -n, x.r = 1/x.r, x.i=1/x.i;
299                 for(u = n; ; ) {
300                         if(u & 01) pow.r *= x.r, pow.i *= x.i;
301                         if(u >>= 1) x.r *= x.r, x.i *= x.i;
302                         else break;
303                 }
304         }
305         _Fcomplex p={pow.r, pow.i};
306         return p;
307 }
308 #else
309 static _Complex float cpow_ui(_Complex float x, integer n) {
310         _Complex float pow=1.0; unsigned long int u;
311         if(n != 0) {
312                 if(n < 0) n = -n, x = 1/x;
313                 for(u = n; ; ) {
314                         if(u & 01) pow *= x;
315                         if(u >>= 1) x *= x;
316                         else break;
317                 }
318         }
319         return pow;
320 }
321 #endif
322 #ifdef _MSC_VER
323 static _Dcomplex zpow_ui(_Dcomplex x, integer n) {
324         _Dcomplex pow={1.0,0.0}; unsigned long int u;
325         if(n != 0) {
326                 if(n < 0) n = -n, x._Val[0] = 1/x._Val[0], x._Val[1] =1/x._Val[1];
327                 for(u = n; ; ) {
328                         if(u & 01) pow._Val[0] *= x._Val[0], pow._Val[1] *= x._Val[1];
329                         if(u >>= 1) x._Val[0] *= x._Val[0], x._Val[1] *= x._Val[1];
330                         else break;
331                 }
332         }
333         _Dcomplex p = {pow._Val[0], pow._Val[1]};
334         return p;
335 }
336 #else
337 static _Complex double zpow_ui(_Complex double x, integer n) {
338         _Complex double pow=1.0; unsigned long int u;
339         if(n != 0) {
340                 if(n < 0) n = -n, x = 1/x;
341                 for(u = n; ; ) {
342                         if(u & 01) pow *= x;
343                         if(u >>= 1) x *= x;
344                         else break;
345                 }
346         }
347         return pow;
348 }
349 #endif
350 static integer pow_ii(integer x, integer n) {
351         integer pow; unsigned long int u;
352         if (n <= 0) {
353                 if (n == 0 || x == 1) pow = 1;
354                 else if (x != -1) pow = x == 0 ? 1/x : 0;
355                 else n = -n;
356         }
357         if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
358                 u = n;
359                 for(pow = 1; ; ) {
360                         if(u & 01) pow *= x;
361                         if(u >>= 1) x *= x;
362                         else break;
363                 }
364         }
365         return pow;
366 }
367 static integer dmaxloc_(double *w, integer s, integer e, integer *n)
368 {
369         double m; integer i, mi;
370         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
371                 if (w[i-1]>m) mi=i ,m=w[i-1];
372         return mi-s+1;
373 }
374 static integer smaxloc_(float *w, integer s, integer e, integer *n)
375 {
376         float m; integer i, mi;
377         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
378                 if (w[i-1]>m) mi=i ,m=w[i-1];
379         return mi-s+1;
380 }
381 static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
382         integer n = *n_, incx = *incx_, incy = *incy_, i;
383 #ifdef _MSC_VER
384         _Fcomplex zdotc = {0.0, 0.0};
385         if (incx == 1 && incy == 1) {
386                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
387                         zdotc._Val[0] += conjf(Cf(&x[i]))._Val[0] * Cf(&y[i])._Val[0];
388                         zdotc._Val[1] += conjf(Cf(&x[i]))._Val[1] * Cf(&y[i])._Val[1];
389                 }
390         } else {
391                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
392                         zdotc._Val[0] += conjf(Cf(&x[i*incx]))._Val[0] * Cf(&y[i*incy])._Val[0];
393                         zdotc._Val[1] += conjf(Cf(&x[i*incx]))._Val[1] * Cf(&y[i*incy])._Val[1];
394                 }
395         }
396         pCf(z) = zdotc;
397 }
398 #else
399         _Complex float zdotc = 0.0;
400         if (incx == 1 && incy == 1) {
401                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
402                         zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
403                 }
404         } else {
405                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
406                         zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
407                 }
408         }
409         pCf(z) = zdotc;
410 }
411 #endif
412 static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
413         integer n = *n_, incx = *incx_, incy = *incy_, i;
414 #ifdef _MSC_VER
415         _Dcomplex zdotc = {0.0, 0.0};
416         if (incx == 1 && incy == 1) {
417                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
418                         zdotc._Val[0] += conj(Cd(&x[i]))._Val[0] * Cd(&y[i])._Val[0];
419                         zdotc._Val[1] += conj(Cd(&x[i]))._Val[1] * Cd(&y[i])._Val[1];
420                 }
421         } else {
422                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
423                         zdotc._Val[0] += conj(Cd(&x[i*incx]))._Val[0] * Cd(&y[i*incy])._Val[0];
424                         zdotc._Val[1] += conj(Cd(&x[i*incx]))._Val[1] * Cd(&y[i*incy])._Val[1];
425                 }
426         }
427         pCd(z) = zdotc;
428 }
429 #else
430         _Complex double zdotc = 0.0;
431         if (incx == 1 && incy == 1) {
432                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
433                         zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
434                 }
435         } else {
436                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
437                         zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
438                 }
439         }
440         pCd(z) = zdotc;
441 }
442 #endif  
443 static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
444         integer n = *n_, incx = *incx_, incy = *incy_, i;
445 #ifdef _MSC_VER
446         _Fcomplex zdotc = {0.0, 0.0};
447         if (incx == 1 && incy == 1) {
448                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
449                         zdotc._Val[0] += Cf(&x[i])._Val[0] * Cf(&y[i])._Val[0];
450                         zdotc._Val[1] += Cf(&x[i])._Val[1] * Cf(&y[i])._Val[1];
451                 }
452         } else {
453                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
454                         zdotc._Val[0] += Cf(&x[i*incx])._Val[0] * Cf(&y[i*incy])._Val[0];
455                         zdotc._Val[1] += Cf(&x[i*incx])._Val[1] * Cf(&y[i*incy])._Val[1];
456                 }
457         }
458         pCf(z) = zdotc;
459 }
460 #else
461         _Complex float zdotc = 0.0;
462         if (incx == 1 && incy == 1) {
463                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
464                         zdotc += Cf(&x[i]) * Cf(&y[i]);
465                 }
466         } else {
467                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
468                         zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
469                 }
470         }
471         pCf(z) = zdotc;
472 }
473 #endif
474 static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
475         integer n = *n_, incx = *incx_, incy = *incy_, i;
476 #ifdef _MSC_VER
477         _Dcomplex zdotc = {0.0, 0.0};
478         if (incx == 1 && incy == 1) {
479                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
480                         zdotc._Val[0] += Cd(&x[i])._Val[0] * Cd(&y[i])._Val[0];
481                         zdotc._Val[1] += Cd(&x[i])._Val[1] * Cd(&y[i])._Val[1];
482                 }
483         } else {
484                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
485                         zdotc._Val[0] += Cd(&x[i*incx])._Val[0] * Cd(&y[i*incy])._Val[0];
486                         zdotc._Val[1] += Cd(&x[i*incx])._Val[1] * Cd(&y[i*incy])._Val[1];
487                 }
488         }
489         pCd(z) = zdotc;
490 }
491 #else
492         _Complex double zdotc = 0.0;
493         if (incx == 1 && incy == 1) {
494                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
495                         zdotc += Cd(&x[i]) * Cd(&y[i]);
496                 }
497         } else {
498                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
499                         zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
500                 }
501         }
502         pCd(z) = zdotc;
503 }
504 #endif
505 /*  -- translated by f2c (version 20000121).
506    You must link the resulting object file with the libraries:
507         -lf2c -lm   (in that order)
508 */
509
510
511
512
513 /* Table of constant values */
514
515 static integer c__1 = 1;
516 static real c_b12 = -1.f;
517 static real c_b14 = 1.f;
518
519 /* > \brief \b SSPRFS */
520
521 /*  =========== DOCUMENTATION =========== */
522
523 /* Online html documentation available at */
524 /*            http://www.netlib.org/lapack/explore-html/ */
525
526 /* > \htmlonly */
527 /* > Download SSPRFS + dependencies */
528 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ssprfs.
529 f"> */
530 /* > [TGZ]</a> */
531 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ssprfs.
532 f"> */
533 /* > [ZIP]</a> */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ssprfs.
535 f"> */
536 /* > [TXT]</a> */
537 /* > \endhtmlonly */
538
539 /*  Definition: */
540 /*  =========== */
541
542 /*       SUBROUTINE SSPRFS( UPLO, N, NRHS, AP, AFP, IPIV, B, LDB, X, LDX, */
543 /*                          FERR, BERR, WORK, IWORK, INFO ) */
544
545 /*       CHARACTER          UPLO */
546 /*       INTEGER            INFO, LDB, LDX, N, NRHS */
547 /*       INTEGER            IPIV( * ), IWORK( * ) */
548 /*       REAL               AFP( * ), AP( * ), B( LDB, * ), BERR( * ), */
549 /*      $                   FERR( * ), WORK( * ), X( LDX, * ) */
550
551
552 /* > \par Purpose: */
553 /*  ============= */
554 /* > */
555 /* > \verbatim */
556 /* > */
557 /* > SSPRFS improves the computed solution to a system of linear */
558 /* > equations when the coefficient matrix is symmetric indefinite */
559 /* > and packed, and provides error bounds and backward error estimates */
560 /* > for the solution. */
561 /* > \endverbatim */
562
563 /*  Arguments: */
564 /*  ========== */
565
566 /* > \param[in] UPLO */
567 /* > \verbatim */
568 /* >          UPLO is CHARACTER*1 */
569 /* >          = 'U':  Upper triangle of A is stored; */
570 /* >          = 'L':  Lower triangle of A is stored. */
571 /* > \endverbatim */
572 /* > */
573 /* > \param[in] N */
574 /* > \verbatim */
575 /* >          N is INTEGER */
576 /* >          The order of the matrix A.  N >= 0. */
577 /* > \endverbatim */
578 /* > */
579 /* > \param[in] NRHS */
580 /* > \verbatim */
581 /* >          NRHS is INTEGER */
582 /* >          The number of right hand sides, i.e., the number of columns */
583 /* >          of the matrices B and X.  NRHS >= 0. */
584 /* > \endverbatim */
585 /* > */
586 /* > \param[in] AP */
587 /* > \verbatim */
588 /* >          AP is REAL array, dimension (N*(N+1)/2) */
589 /* >          The upper or lower triangle of the symmetric matrix A, packed */
590 /* >          columnwise in a linear array.  The j-th column of A is stored */
591 /* >          in the array AP as follows: */
592 /* >          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j; */
593 /* >          if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n. */
594 /* > \endverbatim */
595 /* > */
596 /* > \param[in] AFP */
597 /* > \verbatim */
598 /* >          AFP is REAL array, dimension (N*(N+1)/2) */
599 /* >          The factored form of the matrix A.  AFP contains the block */
600 /* >          diagonal matrix D and the multipliers used to obtain the */
601 /* >          factor U or L from the factorization A = U*D*U**T or */
602 /* >          A = L*D*L**T as computed by SSPTRF, stored as a packed */
603 /* >          triangular matrix. */
604 /* > \endverbatim */
605 /* > */
606 /* > \param[in] IPIV */
607 /* > \verbatim */
608 /* >          IPIV is INTEGER array, dimension (N) */
609 /* >          Details of the interchanges and the block structure of D */
610 /* >          as determined by SSPTRF. */
611 /* > \endverbatim */
612 /* > */
613 /* > \param[in] B */
614 /* > \verbatim */
615 /* >          B is REAL array, dimension (LDB,NRHS) */
616 /* >          The right hand side matrix B. */
617 /* > \endverbatim */
618 /* > */
619 /* > \param[in] LDB */
620 /* > \verbatim */
621 /* >          LDB is INTEGER */
622 /* >          The leading dimension of the array B.  LDB >= f2cmax(1,N). */
623 /* > \endverbatim */
624 /* > */
625 /* > \param[in,out] X */
626 /* > \verbatim */
627 /* >          X is REAL array, dimension (LDX,NRHS) */
628 /* >          On entry, the solution matrix X, as computed by SSPTRS. */
629 /* >          On exit, the improved solution matrix X. */
630 /* > \endverbatim */
631 /* > */
632 /* > \param[in] LDX */
633 /* > \verbatim */
634 /* >          LDX is INTEGER */
635 /* >          The leading dimension of the array X.  LDX >= f2cmax(1,N). */
636 /* > \endverbatim */
637 /* > */
638 /* > \param[out] FERR */
639 /* > \verbatim */
640 /* >          FERR is REAL array, dimension (NRHS) */
641 /* >          The estimated forward error bound for each solution vector */
642 /* >          X(j) (the j-th column of the solution matrix X). */
643 /* >          If XTRUE is the true solution corresponding to X(j), FERR(j) */
644 /* >          is an estimated upper bound for the magnitude of the largest */
645 /* >          element in (X(j) - XTRUE) divided by the magnitude of the */
646 /* >          largest element in X(j).  The estimate is as reliable as */
647 /* >          the estimate for RCOND, and is almost always a slight */
648 /* >          overestimate of the true error. */
649 /* > \endverbatim */
650 /* > */
651 /* > \param[out] BERR */
652 /* > \verbatim */
653 /* >          BERR is REAL array, dimension (NRHS) */
654 /* >          The componentwise relative backward error of each solution */
655 /* >          vector X(j) (i.e., the smallest relative change in */
656 /* >          any element of A or B that makes X(j) an exact solution). */
657 /* > \endverbatim */
658 /* > */
659 /* > \param[out] WORK */
660 /* > \verbatim */
661 /* >          WORK is REAL array, dimension (3*N) */
662 /* > \endverbatim */
663 /* > */
664 /* > \param[out] IWORK */
665 /* > \verbatim */
666 /* >          IWORK is INTEGER array, dimension (N) */
667 /* > \endverbatim */
668 /* > */
669 /* > \param[out] INFO */
670 /* > \verbatim */
671 /* >          INFO is INTEGER */
672 /* >          = 0:  successful exit */
673 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value */
674 /* > \endverbatim */
675
676 /* > \par Internal Parameters: */
677 /*  ========================= */
678 /* > */
679 /* > \verbatim */
680 /* >  ITMAX is the maximum number of steps of iterative refinement. */
681 /* > \endverbatim */
682
683 /*  Authors: */
684 /*  ======== */
685
686 /* > \author Univ. of Tennessee */
687 /* > \author Univ. of California Berkeley */
688 /* > \author Univ. of Colorado Denver */
689 /* > \author NAG Ltd. */
690
691 /* > \date December 2016 */
692
693 /* > \ingroup realOTHERcomputational */
694
695 /*  ===================================================================== */
696 /* Subroutine */ int ssprfs_(char *uplo, integer *n, integer *nrhs, real *ap, 
697         real *afp, integer *ipiv, real *b, integer *ldb, real *x, integer *
698         ldx, real *ferr, real *berr, real *work, integer *iwork, integer *
699         info)
700 {
701     /* System generated locals */
702     integer b_dim1, b_offset, x_dim1, x_offset, i__1, i__2, i__3;
703     real r__1, r__2, r__3;
704
705     /* Local variables */
706     integer kase;
707     real safe1, safe2;
708     integer i__, j, k;
709     real s;
710     extern logical lsame_(char *, char *);
711     integer isave[3], count;
712     logical upper;
713     extern /* Subroutine */ int scopy_(integer *, real *, integer *, real *, 
714             integer *), saxpy_(integer *, real *, real *, integer *, real *, 
715             integer *), sspmv_(char *, integer *, real *, real *, real *, 
716             integer *, real *, real *, integer *), slacn2_(integer *, 
717             real *, real *, integer *, real *, integer *, integer *);
718     integer ik, kk;
719     real xk;
720     extern real slamch_(char *);
721     integer nz;
722     real safmin;
723     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
724     real lstres;
725     extern /* Subroutine */ int ssptrs_(char *, integer *, integer *, real *, 
726             integer *, real *, integer *, integer *);
727     real eps;
728
729
730 /*  -- LAPACK computational routine (version 3.7.0) -- */
731 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
732 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
733 /*     December 2016 */
734
735
736 /*  ===================================================================== */
737
738
739 /*     Test the input parameters. */
740
741     /* Parameter adjustments */
742     --ap;
743     --afp;
744     --ipiv;
745     b_dim1 = *ldb;
746     b_offset = 1 + b_dim1 * 1;
747     b -= b_offset;
748     x_dim1 = *ldx;
749     x_offset = 1 + x_dim1 * 1;
750     x -= x_offset;
751     --ferr;
752     --berr;
753     --work;
754     --iwork;
755
756     /* Function Body */
757     *info = 0;
758     upper = lsame_(uplo, "U");
759     if (! upper && ! lsame_(uplo, "L")) {
760         *info = -1;
761     } else if (*n < 0) {
762         *info = -2;
763     } else if (*nrhs < 0) {
764         *info = -3;
765     } else if (*ldb < f2cmax(1,*n)) {
766         *info = -8;
767     } else if (*ldx < f2cmax(1,*n)) {
768         *info = -10;
769     }
770     if (*info != 0) {
771         i__1 = -(*info);
772         xerbla_("SSPRFS", &i__1, (ftnlen)6);
773         return 0;
774     }
775
776 /*     Quick return if possible */
777
778     if (*n == 0 || *nrhs == 0) {
779         i__1 = *nrhs;
780         for (j = 1; j <= i__1; ++j) {
781             ferr[j] = 0.f;
782             berr[j] = 0.f;
783 /* L10: */
784         }
785         return 0;
786     }
787
788 /*     NZ = maximum number of nonzero elements in each row of A, plus 1 */
789
790     nz = *n + 1;
791     eps = slamch_("Epsilon");
792     safmin = slamch_("Safe minimum");
793     safe1 = nz * safmin;
794     safe2 = safe1 / eps;
795
796 /*     Do for each right hand side */
797
798     i__1 = *nrhs;
799     for (j = 1; j <= i__1; ++j) {
800
801         count = 1;
802         lstres = 3.f;
803 L20:
804
805 /*        Loop until stopping criterion is satisfied. */
806
807 /*        Compute residual R = B - A * X */
808
809         scopy_(n, &b[j * b_dim1 + 1], &c__1, &work[*n + 1], &c__1);
810         sspmv_(uplo, n, &c_b12, &ap[1], &x[j * x_dim1 + 1], &c__1, &c_b14, &
811                 work[*n + 1], &c__1);
812
813 /*        Compute componentwise relative backward error from formula */
814
815 /*        f2cmax(i) ( abs(R(i)) / ( abs(A)*abs(X) + abs(B) )(i) ) */
816
817 /*        where abs(Z) is the componentwise absolute value of the matrix */
818 /*        or vector Z.  If the i-th component of the denominator is less */
819 /*        than SAFE2, then SAFE1 is added to the i-th components of the */
820 /*        numerator and denominator before dividing. */
821
822         i__2 = *n;
823         for (i__ = 1; i__ <= i__2; ++i__) {
824             work[i__] = (r__1 = b[i__ + j * b_dim1], abs(r__1));
825 /* L30: */
826         }
827
828 /*        Compute abs(A)*abs(X) + abs(B). */
829
830         kk = 1;
831         if (upper) {
832             i__2 = *n;
833             for (k = 1; k <= i__2; ++k) {
834                 s = 0.f;
835                 xk = (r__1 = x[k + j * x_dim1], abs(r__1));
836                 ik = kk;
837                 i__3 = k - 1;
838                 for (i__ = 1; i__ <= i__3; ++i__) {
839                     work[i__] += (r__1 = ap[ik], abs(r__1)) * xk;
840                     s += (r__1 = ap[ik], abs(r__1)) * (r__2 = x[i__ + j * 
841                             x_dim1], abs(r__2));
842                     ++ik;
843 /* L40: */
844                 }
845                 work[k] = work[k] + (r__1 = ap[kk + k - 1], abs(r__1)) * xk + 
846                         s;
847                 kk += k;
848 /* L50: */
849             }
850         } else {
851             i__2 = *n;
852             for (k = 1; k <= i__2; ++k) {
853                 s = 0.f;
854                 xk = (r__1 = x[k + j * x_dim1], abs(r__1));
855                 work[k] += (r__1 = ap[kk], abs(r__1)) * xk;
856                 ik = kk + 1;
857                 i__3 = *n;
858                 for (i__ = k + 1; i__ <= i__3; ++i__) {
859                     work[i__] += (r__1 = ap[ik], abs(r__1)) * xk;
860                     s += (r__1 = ap[ik], abs(r__1)) * (r__2 = x[i__ + j * 
861                             x_dim1], abs(r__2));
862                     ++ik;
863 /* L60: */
864                 }
865                 work[k] += s;
866                 kk += *n - k + 1;
867 /* L70: */
868             }
869         }
870         s = 0.f;
871         i__2 = *n;
872         for (i__ = 1; i__ <= i__2; ++i__) {
873             if (work[i__] > safe2) {
874 /* Computing MAX */
875                 r__2 = s, r__3 = (r__1 = work[*n + i__], abs(r__1)) / work[
876                         i__];
877                 s = f2cmax(r__2,r__3);
878             } else {
879 /* Computing MAX */
880                 r__2 = s, r__3 = ((r__1 = work[*n + i__], abs(r__1)) + safe1) 
881                         / (work[i__] + safe1);
882                 s = f2cmax(r__2,r__3);
883             }
884 /* L80: */
885         }
886         berr[j] = s;
887
888 /*        Test stopping criterion. Continue iterating if */
889 /*           1) The residual BERR(J) is larger than machine epsilon, and */
890 /*           2) BERR(J) decreased by at least a factor of 2 during the */
891 /*              last iteration, and */
892 /*           3) At most ITMAX iterations tried. */
893
894         if (berr[j] > eps && berr[j] * 2.f <= lstres && count <= 5) {
895
896 /*           Update solution and try again. */
897
898             ssptrs_(uplo, n, &c__1, &afp[1], &ipiv[1], &work[*n + 1], n, info);
899             saxpy_(n, &c_b14, &work[*n + 1], &c__1, &x[j * x_dim1 + 1], &c__1)
900                     ;
901             lstres = berr[j];
902             ++count;
903             goto L20;
904         }
905
906 /*        Bound error from formula */
907
908 /*        norm(X - XTRUE) / norm(X) .le. FERR = */
909 /*        norm( abs(inv(A))* */
910 /*           ( abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) ))) / norm(X) */
911
912 /*        where */
913 /*          norm(Z) is the magnitude of the largest component of Z */
914 /*          inv(A) is the inverse of A */
915 /*          abs(Z) is the componentwise absolute value of the matrix or */
916 /*             vector Z */
917 /*          NZ is the maximum number of nonzeros in any row of A, plus 1 */
918 /*          EPS is machine epsilon */
919
920 /*        The i-th component of abs(R)+NZ*EPS*(abs(A)*abs(X)+abs(B)) */
921 /*        is incremented by SAFE1 if the i-th component of */
922 /*        abs(A)*abs(X) + abs(B) is less than SAFE2. */
923
924 /*        Use SLACN2 to estimate the infinity-norm of the matrix */
925 /*           inv(A) * diag(W), */
926 /*        where W = abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) ))) */
927
928         i__2 = *n;
929         for (i__ = 1; i__ <= i__2; ++i__) {
930             if (work[i__] > safe2) {
931                 work[i__] = (r__1 = work[*n + i__], abs(r__1)) + nz * eps * 
932                         work[i__];
933             } else {
934                 work[i__] = (r__1 = work[*n + i__], abs(r__1)) + nz * eps * 
935                         work[i__] + safe1;
936             }
937 /* L90: */
938         }
939
940         kase = 0;
941 L100:
942         slacn2_(n, &work[(*n << 1) + 1], &work[*n + 1], &iwork[1], &ferr[j], &
943                 kase, isave);
944         if (kase != 0) {
945             if (kase == 1) {
946
947 /*              Multiply by diag(W)*inv(A**T). */
948
949                 ssptrs_(uplo, n, &c__1, &afp[1], &ipiv[1], &work[*n + 1], n, 
950                         info);
951                 i__2 = *n;
952                 for (i__ = 1; i__ <= i__2; ++i__) {
953                     work[*n + i__] = work[i__] * work[*n + i__];
954 /* L110: */
955                 }
956             } else if (kase == 2) {
957
958 /*              Multiply by inv(A)*diag(W). */
959
960                 i__2 = *n;
961                 for (i__ = 1; i__ <= i__2; ++i__) {
962                     work[*n + i__] = work[i__] * work[*n + i__];
963 /* L120: */
964                 }
965                 ssptrs_(uplo, n, &c__1, &afp[1], &ipiv[1], &work[*n + 1], n, 
966                         info);
967             }
968             goto L100;
969         }
970
971 /*        Normalize error. */
972
973         lstres = 0.f;
974         i__2 = *n;
975         for (i__ = 1; i__ <= i__2; ++i__) {
976 /* Computing MAX */
977             r__2 = lstres, r__3 = (r__1 = x[i__ + j * x_dim1], abs(r__1));
978             lstres = f2cmax(r__2,r__3);
979 /* L130: */
980         }
981         if (lstres != 0.f) {
982             ferr[j] /= lstres;
983         }
984
985 /* L140: */
986     }
987
988     return 0;
989
990 /*     End of SSPRFS */
991
992 } /* ssprfs_ */
993