C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / sgbrfs.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_b15 = -1.f;
517 static real c_b17 = 1.f;
518
519 /* > \brief \b SGBRFS */
520
521 /*  =========== DOCUMENTATION =========== */
522
523 /* Online html documentation available at */
524 /*            http://www.netlib.org/lapack/explore-html/ */
525
526 /* > \htmlonly */
527 /* > Download SGBRFS + dependencies */
528 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgbrfs.
529 f"> */
530 /* > [TGZ]</a> */
531 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgbrfs.
532 f"> */
533 /* > [ZIP]</a> */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgbrfs.
535 f"> */
536 /* > [TXT]</a> */
537 /* > \endhtmlonly */
538
539 /*  Definition: */
540 /*  =========== */
541
542 /*       SUBROUTINE SGBRFS( TRANS, N, KL, KU, NRHS, AB, LDAB, AFB, LDAFB, */
543 /*                          IPIV, B, LDB, X, LDX, FERR, BERR, WORK, IWORK, */
544 /*                          INFO ) */
545
546 /*       CHARACTER          TRANS */
547 /*       INTEGER            INFO, KL, KU, LDAB, LDAFB, LDB, LDX, N, NRHS */
548 /*       INTEGER            IPIV( * ), IWORK( * ) */
549 /*       REAL               AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ), */
550 /*      $                   BERR( * ), FERR( * ), WORK( * ), X( LDX, * ) */
551
552
553 /* > \par Purpose: */
554 /*  ============= */
555 /* > */
556 /* > \verbatim */
557 /* > */
558 /* > SGBRFS improves the computed solution to a system of linear */
559 /* > equations when the coefficient matrix is banded, and provides */
560 /* > error bounds and backward error estimates for the solution. */
561 /* > \endverbatim */
562
563 /*  Arguments: */
564 /*  ========== */
565
566 /* > \param[in] TRANS */
567 /* > \verbatim */
568 /* >          TRANS is CHARACTER*1 */
569 /* >          Specifies the form of the system of equations: */
570 /* >          = 'N':  A * X = B     (No transpose) */
571 /* >          = 'T':  A**T * X = B  (Transpose) */
572 /* >          = 'C':  A**H * X = B  (Conjugate transpose = Transpose) */
573 /* > \endverbatim */
574 /* > */
575 /* > \param[in] N */
576 /* > \verbatim */
577 /* >          N is INTEGER */
578 /* >          The order of the matrix A.  N >= 0. */
579 /* > \endverbatim */
580 /* > */
581 /* > \param[in] KL */
582 /* > \verbatim */
583 /* >          KL is INTEGER */
584 /* >          The number of subdiagonals within the band of A.  KL >= 0. */
585 /* > \endverbatim */
586 /* > */
587 /* > \param[in] KU */
588 /* > \verbatim */
589 /* >          KU is INTEGER */
590 /* >          The number of superdiagonals within the band of A.  KU >= 0. */
591 /* > \endverbatim */
592 /* > */
593 /* > \param[in] NRHS */
594 /* > \verbatim */
595 /* >          NRHS is INTEGER */
596 /* >          The number of right hand sides, i.e., the number of columns */
597 /* >          of the matrices B and X.  NRHS >= 0. */
598 /* > \endverbatim */
599 /* > */
600 /* > \param[in] AB */
601 /* > \verbatim */
602 /* >          AB is REAL array, dimension (LDAB,N) */
603 /* >          The original band matrix A, stored in rows 1 to KL+KU+1. */
604 /* >          The j-th column of A is stored in the j-th column of the */
605 /* >          array AB as follows: */
606 /* >          AB(ku+1+i-j,j) = A(i,j) for f2cmax(1,j-ku)<=i<=f2cmin(n,j+kl). */
607 /* > \endverbatim */
608 /* > */
609 /* > \param[in] LDAB */
610 /* > \verbatim */
611 /* >          LDAB is INTEGER */
612 /* >          The leading dimension of the array AB.  LDAB >= KL+KU+1. */
613 /* > \endverbatim */
614 /* > */
615 /* > \param[in] AFB */
616 /* > \verbatim */
617 /* >          AFB is REAL array, dimension (LDAFB,N) */
618 /* >          Details of the LU factorization of the band matrix A, as */
619 /* >          computed by SGBTRF.  U is stored as an upper triangular band */
620 /* >          matrix with KL+KU superdiagonals in rows 1 to KL+KU+1, and */
621 /* >          the multipliers used during the factorization are stored in */
622 /* >          rows KL+KU+2 to 2*KL+KU+1. */
623 /* > \endverbatim */
624 /* > */
625 /* > \param[in] LDAFB */
626 /* > \verbatim */
627 /* >          LDAFB is INTEGER */
628 /* >          The leading dimension of the array AFB.  LDAFB >= 2*KL*KU+1. */
629 /* > \endverbatim */
630 /* > */
631 /* > \param[in] IPIV */
632 /* > \verbatim */
633 /* >          IPIV is INTEGER array, dimension (N) */
634 /* >          The pivot indices from SGBTRF; for 1<=i<=N, row i of the */
635 /* >          matrix was interchanged with row IPIV(i). */
636 /* > \endverbatim */
637 /* > */
638 /* > \param[in] B */
639 /* > \verbatim */
640 /* >          B is REAL array, dimension (LDB,NRHS) */
641 /* >          The right hand side matrix B. */
642 /* > \endverbatim */
643 /* > */
644 /* > \param[in] LDB */
645 /* > \verbatim */
646 /* >          LDB is INTEGER */
647 /* >          The leading dimension of the array B.  LDB >= f2cmax(1,N). */
648 /* > \endverbatim */
649 /* > */
650 /* > \param[in,out] X */
651 /* > \verbatim */
652 /* >          X is REAL array, dimension (LDX,NRHS) */
653 /* >          On entry, the solution matrix X, as computed by SGBTRS. */
654 /* >          On exit, the improved solution matrix X. */
655 /* > \endverbatim */
656 /* > */
657 /* > \param[in] LDX */
658 /* > \verbatim */
659 /* >          LDX is INTEGER */
660 /* >          The leading dimension of the array X.  LDX >= f2cmax(1,N). */
661 /* > \endverbatim */
662 /* > */
663 /* > \param[out] FERR */
664 /* > \verbatim */
665 /* >          FERR is REAL array, dimension (NRHS) */
666 /* >          The estimated forward error bound for each solution vector */
667 /* >          X(j) (the j-th column of the solution matrix X). */
668 /* >          If XTRUE is the true solution corresponding to X(j), FERR(j) */
669 /* >          is an estimated upper bound for the magnitude of the largest */
670 /* >          element in (X(j) - XTRUE) divided by the magnitude of the */
671 /* >          largest element in X(j).  The estimate is as reliable as */
672 /* >          the estimate for RCOND, and is almost always a slight */
673 /* >          overestimate of the true error. */
674 /* > \endverbatim */
675 /* > */
676 /* > \param[out] BERR */
677 /* > \verbatim */
678 /* >          BERR is REAL array, dimension (NRHS) */
679 /* >          The componentwise relative backward error of each solution */
680 /* >          vector X(j) (i.e., the smallest relative change in */
681 /* >          any element of A or B that makes X(j) an exact solution). */
682 /* > \endverbatim */
683 /* > */
684 /* > \param[out] WORK */
685 /* > \verbatim */
686 /* >          WORK is REAL array, dimension (3*N) */
687 /* > \endverbatim */
688 /* > */
689 /* > \param[out] IWORK */
690 /* > \verbatim */
691 /* >          IWORK is INTEGER array, dimension (N) */
692 /* > \endverbatim */
693 /* > */
694 /* > \param[out] INFO */
695 /* > \verbatim */
696 /* >          INFO is INTEGER */
697 /* >          = 0:  successful exit */
698 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value */
699 /* > \endverbatim */
700
701 /* > \par Internal Parameters: */
702 /*  ========================= */
703 /* > */
704 /* > \verbatim */
705 /* >  ITMAX is the maximum number of steps of iterative refinement. */
706 /* > \endverbatim */
707
708 /*  Authors: */
709 /*  ======== */
710
711 /* > \author Univ. of Tennessee */
712 /* > \author Univ. of California Berkeley */
713 /* > \author Univ. of Colorado Denver */
714 /* > \author NAG Ltd. */
715
716 /* > \date December 2016 */
717
718 /* > \ingroup realGBcomputational */
719
720 /*  ===================================================================== */
721 /* Subroutine */ int sgbrfs_(char *trans, integer *n, integer *kl, integer *
722         ku, integer *nrhs, real *ab, integer *ldab, real *afb, integer *ldafb,
723          integer *ipiv, real *b, integer *ldb, real *x, integer *ldx, real *
724         ferr, real *berr, real *work, integer *iwork, integer *info)
725 {
726     /* System generated locals */
727     integer ab_dim1, ab_offset, afb_dim1, afb_offset, b_dim1, b_offset, 
728             x_dim1, x_offset, i__1, i__2, i__3, i__4, i__5, i__6, i__7;
729     real r__1, r__2, r__3;
730
731     /* Local variables */
732     integer kase;
733     real safe1, safe2;
734     integer i__, j, k;
735     real s;
736     extern logical lsame_(char *, char *);
737     integer isave[3];
738     extern /* Subroutine */ int sgbmv_(char *, integer *, integer *, integer *
739             , integer *, real *, real *, integer *, real *, integer *, real *,
740              real *, integer *);
741     integer count;
742     extern /* Subroutine */ int scopy_(integer *, real *, integer *, real *, 
743             integer *), saxpy_(integer *, real *, real *, integer *, real *, 
744             integer *), slacn2_(integer *, real *, real *, integer *, real *, 
745             integer *, integer *);
746     integer kk;
747     real xk;
748     extern real slamch_(char *);
749     integer nz;
750     real safmin;
751     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
752     logical notran;
753     extern /* Subroutine */ int sgbtrs_(char *, integer *, integer *, integer 
754             *, integer *, real *, integer *, integer *, real *, integer *, 
755             integer *);
756     char transt[1];
757     real lstres, eps;
758
759
760 /*  -- LAPACK computational routine (version 3.7.0) -- */
761 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
762 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
763 /*     December 2016 */
764
765
766 /*  ===================================================================== */
767
768
769 /*     Test the input parameters. */
770
771     /* Parameter adjustments */
772     ab_dim1 = *ldab;
773     ab_offset = 1 + ab_dim1 * 1;
774     ab -= ab_offset;
775     afb_dim1 = *ldafb;
776     afb_offset = 1 + afb_dim1 * 1;
777     afb -= afb_offset;
778     --ipiv;
779     b_dim1 = *ldb;
780     b_offset = 1 + b_dim1 * 1;
781     b -= b_offset;
782     x_dim1 = *ldx;
783     x_offset = 1 + x_dim1 * 1;
784     x -= x_offset;
785     --ferr;
786     --berr;
787     --work;
788     --iwork;
789
790     /* Function Body */
791     *info = 0;
792     notran = lsame_(trans, "N");
793     if (! notran && ! lsame_(trans, "T") && ! lsame_(
794             trans, "C")) {
795         *info = -1;
796     } else if (*n < 0) {
797         *info = -2;
798     } else if (*kl < 0) {
799         *info = -3;
800     } else if (*ku < 0) {
801         *info = -4;
802     } else if (*nrhs < 0) {
803         *info = -5;
804     } else if (*ldab < *kl + *ku + 1) {
805         *info = -7;
806     } else if (*ldafb < (*kl << 1) + *ku + 1) {
807         *info = -9;
808     } else if (*ldb < f2cmax(1,*n)) {
809         *info = -12;
810     } else if (*ldx < f2cmax(1,*n)) {
811         *info = -14;
812     }
813     if (*info != 0) {
814         i__1 = -(*info);
815         xerbla_("SGBRFS", &i__1, (ftnlen)6);
816         return 0;
817     }
818
819 /*     Quick return if possible */
820
821     if (*n == 0 || *nrhs == 0) {
822         i__1 = *nrhs;
823         for (j = 1; j <= i__1; ++j) {
824             ferr[j] = 0.f;
825             berr[j] = 0.f;
826 /* L10: */
827         }
828         return 0;
829     }
830
831     if (notran) {
832         *(unsigned char *)transt = 'T';
833     } else {
834         *(unsigned char *)transt = 'N';
835     }
836
837 /*     NZ = maximum number of nonzero elements in each row of A, plus 1 */
838
839 /* Computing MIN */
840     i__1 = *kl + *ku + 2, i__2 = *n + 1;
841     nz = f2cmin(i__1,i__2);
842     eps = slamch_("Epsilon");
843     safmin = slamch_("Safe minimum");
844     safe1 = nz * safmin;
845     safe2 = safe1 / eps;
846
847 /*     Do for each right hand side */
848
849     i__1 = *nrhs;
850     for (j = 1; j <= i__1; ++j) {
851
852         count = 1;
853         lstres = 3.f;
854 L20:
855
856 /*        Loop until stopping criterion is satisfied. */
857
858 /*        Compute residual R = B - op(A) * X, */
859 /*        where op(A) = A, A**T, or A**H, depending on TRANS. */
860
861         scopy_(n, &b[j * b_dim1 + 1], &c__1, &work[*n + 1], &c__1);
862         sgbmv_(trans, n, n, kl, ku, &c_b15, &ab[ab_offset], ldab, &x[j * 
863                 x_dim1 + 1], &c__1, &c_b17, &work[*n + 1], &c__1);
864
865 /*        Compute componentwise relative backward error from formula */
866
867 /*        f2cmax(i) ( abs(R(i)) / ( abs(op(A))*abs(X) + abs(B) )(i) ) */
868
869 /*        where abs(Z) is the componentwise absolute value of the matrix */
870 /*        or vector Z.  If the i-th component of the denominator is less */
871 /*        than SAFE2, then SAFE1 is added to the i-th components of the */
872 /*        numerator and denominator before dividing. */
873
874         i__2 = *n;
875         for (i__ = 1; i__ <= i__2; ++i__) {
876             work[i__] = (r__1 = b[i__ + j * b_dim1], abs(r__1));
877 /* L30: */
878         }
879
880 /*        Compute abs(op(A))*abs(X) + abs(B). */
881
882         if (notran) {
883             i__2 = *n;
884             for (k = 1; k <= i__2; ++k) {
885                 kk = *ku + 1 - k;
886                 xk = (r__1 = x[k + j * x_dim1], abs(r__1));
887 /* Computing MAX */
888                 i__3 = 1, i__4 = k - *ku;
889 /* Computing MIN */
890                 i__6 = *n, i__7 = k + *kl;
891                 i__5 = f2cmin(i__6,i__7);
892                 for (i__ = f2cmax(i__3,i__4); i__ <= i__5; ++i__) {
893                     work[i__] += (r__1 = ab[kk + i__ + k * ab_dim1], abs(r__1)
894                             ) * xk;
895 /* L40: */
896                 }
897 /* L50: */
898             }
899         } else {
900             i__2 = *n;
901             for (k = 1; k <= i__2; ++k) {
902                 s = 0.f;
903                 kk = *ku + 1 - k;
904 /* Computing MAX */
905                 i__5 = 1, i__3 = k - *ku;
906 /* Computing MIN */
907                 i__6 = *n, i__7 = k + *kl;
908                 i__4 = f2cmin(i__6,i__7);
909                 for (i__ = f2cmax(i__5,i__3); i__ <= i__4; ++i__) {
910                     s += (r__1 = ab[kk + i__ + k * ab_dim1], abs(r__1)) * (
911                             r__2 = x[i__ + j * x_dim1], abs(r__2));
912 /* L60: */
913                 }
914                 work[k] += s;
915 /* L70: */
916             }
917         }
918         s = 0.f;
919         i__2 = *n;
920         for (i__ = 1; i__ <= i__2; ++i__) {
921             if (work[i__] > safe2) {
922 /* Computing MAX */
923                 r__2 = s, r__3 = (r__1 = work[*n + i__], abs(r__1)) / work[
924                         i__];
925                 s = f2cmax(r__2,r__3);
926             } else {
927 /* Computing MAX */
928                 r__2 = s, r__3 = ((r__1 = work[*n + i__], abs(r__1)) + safe1) 
929                         / (work[i__] + safe1);
930                 s = f2cmax(r__2,r__3);
931             }
932 /* L80: */
933         }
934         berr[j] = s;
935
936 /*        Test stopping criterion. Continue iterating if */
937 /*           1) The residual BERR(J) is larger than machine epsilon, and */
938 /*           2) BERR(J) decreased by at least a factor of 2 during the */
939 /*              last iteration, and */
940 /*           3) At most ITMAX iterations tried. */
941
942         if (berr[j] > eps && berr[j] * 2.f <= lstres && count <= 5) {
943
944 /*           Update solution and try again. */
945
946             sgbtrs_(trans, n, kl, ku, &c__1, &afb[afb_offset], ldafb, &ipiv[1]
947                     , &work[*n + 1], n, info);
948             saxpy_(n, &c_b17, &work[*n + 1], &c__1, &x[j * x_dim1 + 1], &c__1)
949                     ;
950             lstres = berr[j];
951             ++count;
952             goto L20;
953         }
954
955 /*        Bound error from formula */
956
957 /*        norm(X - XTRUE) / norm(X) .le. FERR = */
958 /*        norm( abs(inv(op(A)))* */
959 /*           ( abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) ))) / norm(X) */
960
961 /*        where */
962 /*          norm(Z) is the magnitude of the largest component of Z */
963 /*          inv(op(A)) is the inverse of op(A) */
964 /*          abs(Z) is the componentwise absolute value of the matrix or */
965 /*             vector Z */
966 /*          NZ is the maximum number of nonzeros in any row of A, plus 1 */
967 /*          EPS is machine epsilon */
968
969 /*        The i-th component of abs(R)+NZ*EPS*(abs(op(A))*abs(X)+abs(B)) */
970 /*        is incremented by SAFE1 if the i-th component of */
971 /*        abs(op(A))*abs(X) + abs(B) is less than SAFE2. */
972
973 /*        Use SLACN2 to estimate the infinity-norm of the matrix */
974 /*           inv(op(A)) * diag(W), */
975 /*        where W = abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) ))) */
976
977         i__2 = *n;
978         for (i__ = 1; i__ <= i__2; ++i__) {
979             if (work[i__] > safe2) {
980                 work[i__] = (r__1 = work[*n + i__], abs(r__1)) + nz * eps * 
981                         work[i__];
982             } else {
983                 work[i__] = (r__1 = work[*n + i__], abs(r__1)) + nz * eps * 
984                         work[i__] + safe1;
985             }
986 /* L90: */
987         }
988
989         kase = 0;
990 L100:
991         slacn2_(n, &work[(*n << 1) + 1], &work[*n + 1], &iwork[1], &ferr[j], &
992                 kase, isave);
993         if (kase != 0) {
994             if (kase == 1) {
995
996 /*              Multiply by diag(W)*inv(op(A)**T). */
997
998                 sgbtrs_(transt, n, kl, ku, &c__1, &afb[afb_offset], ldafb, &
999                         ipiv[1], &work[*n + 1], n, info);
1000                 i__2 = *n;
1001                 for (i__ = 1; i__ <= i__2; ++i__) {
1002                     work[*n + i__] *= work[i__];
1003 /* L110: */
1004                 }
1005             } else {
1006
1007 /*              Multiply by inv(op(A))*diag(W). */
1008
1009                 i__2 = *n;
1010                 for (i__ = 1; i__ <= i__2; ++i__) {
1011                     work[*n + i__] *= work[i__];
1012 /* L120: */
1013                 }
1014                 sgbtrs_(trans, n, kl, ku, &c__1, &afb[afb_offset], ldafb, &
1015                         ipiv[1], &work[*n + 1], n, info);
1016             }
1017             goto L100;
1018         }
1019
1020 /*        Normalize error. */
1021
1022         lstres = 0.f;
1023         i__2 = *n;
1024         for (i__ = 1; i__ <= i__2; ++i__) {
1025 /* Computing MAX */
1026             r__2 = lstres, r__3 = (r__1 = x[i__ + j * x_dim1], abs(r__1));
1027             lstres = f2cmax(r__2,r__3);
1028 /* L130: */
1029         }
1030         if (lstres != 0.f) {
1031             ferr[j] /= lstres;
1032         }
1033
1034 /* L140: */
1035     }
1036
1037     return 0;
1038
1039 /*     End of SGBRFS */
1040
1041 } /* sgbrfs_ */
1042