C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / cla_gbrcond_c.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
517 /* > \brief \b CLA_GBRCOND_C computes the infinity norm condition number of op(A)*inv(diag(c)) for general ban
518 ded matrices. */
519
520 /*  =========== DOCUMENTATION =========== */
521
522 /* Online html documentation available at */
523 /*            http://www.netlib.org/lapack/explore-html/ */
524
525 /* > \htmlonly */
526 /* > Download CLA_GBRCOND_C + dependencies */
527 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cla_gbr
528 cond_c.f"> */
529 /* > [TGZ]</a> */
530 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cla_gbr
531 cond_c.f"> */
532 /* > [ZIP]</a> */
533 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cla_gbr
534 cond_c.f"> */
535 /* > [TXT]</a> */
536 /* > \endhtmlonly */
537
538 /*  Definition: */
539 /*  =========== */
540
541 /*       REAL FUNCTION CLA_GBRCOND_C( TRANS, N, KL, KU, AB, LDAB, AFB, */
542 /*                                    LDAFB, IPIV, C, CAPPLY, INFO, WORK, */
543 /*                                    RWORK ) */
544
545 /*       CHARACTER          TRANS */
546 /*       LOGICAL            CAPPLY */
547 /*       INTEGER            N, KL, KU, KD, KE, LDAB, LDAFB, INFO */
548 /*       INTEGER            IPIV( * ) */
549 /*       COMPLEX            AB( LDAB, * ), AFB( LDAFB, * ), WORK( * ) */
550 /*       REAL               C( * ), RWORK( * ) */
551
552
553 /* > \par Purpose: */
554 /*  ============= */
555 /* > */
556 /* > \verbatim */
557 /* > */
558 /* >    CLA_GBRCOND_C Computes the infinity norm condition number of */
559 /* >    op(A) * inv(diag(C)) where C is a REAL vector. */
560 /* > \endverbatim */
561
562 /*  Arguments: */
563 /*  ========== */
564
565 /* > \param[in] TRANS */
566 /* > \verbatim */
567 /* >          TRANS is CHARACTER*1 */
568 /* >     Specifies the form of the system of equations: */
569 /* >       = 'N':  A * X = B     (No transpose) */
570 /* >       = 'T':  A**T * X = B  (Transpose) */
571 /* >       = 'C':  A**H * X = B  (Conjugate Transpose = Transpose) */
572 /* > \endverbatim */
573 /* > */
574 /* > \param[in] N */
575 /* > \verbatim */
576 /* >          N is INTEGER */
577 /* >     The number of linear equations, i.e., the order of the */
578 /* >     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] AB */
594 /* > \verbatim */
595 /* >          AB is COMPLEX array, dimension (LDAB,N) */
596 /* >     On entry, the matrix A in band storage, in rows 1 to KL+KU+1. */
597 /* >     The j-th column of A is stored in the j-th column of the */
598 /* >     array AB as follows: */
599 /* >     AB(KU+1+i-j,j) = A(i,j) for f2cmax(1,j-KU)<=i<=f2cmin(N,j+kl) */
600 /* > \endverbatim */
601 /* > */
602 /* > \param[in] LDAB */
603 /* > \verbatim */
604 /* >          LDAB is INTEGER */
605 /* >     The leading dimension of the array AB.  LDAB >= KL+KU+1. */
606 /* > \endverbatim */
607 /* > */
608 /* > \param[in] AFB */
609 /* > \verbatim */
610 /* >          AFB is COMPLEX array, dimension (LDAFB,N) */
611 /* >     Details of the LU factorization of the band matrix A, as */
612 /* >     computed by CGBTRF.  U is stored as an upper triangular */
613 /* >     band matrix with KL+KU superdiagonals in rows 1 to KL+KU+1, */
614 /* >     and the multipliers used during the factorization are stored */
615 /* >     in rows KL+KU+2 to 2*KL+KU+1. */
616 /* > \endverbatim */
617 /* > */
618 /* > \param[in] LDAFB */
619 /* > \verbatim */
620 /* >          LDAFB is INTEGER */
621 /* >     The leading dimension of the array AFB.  LDAFB >= 2*KL+KU+1. */
622 /* > \endverbatim */
623 /* > */
624 /* > \param[in] IPIV */
625 /* > \verbatim */
626 /* >          IPIV is INTEGER array, dimension (N) */
627 /* >     The pivot indices from the factorization A = P*L*U */
628 /* >     as computed by CGBTRF; row i of the matrix was interchanged */
629 /* >     with row IPIV(i). */
630 /* > \endverbatim */
631 /* > */
632 /* > \param[in] C */
633 /* > \verbatim */
634 /* >          C is REAL array, dimension (N) */
635 /* >     The vector C in the formula op(A) * inv(diag(C)). */
636 /* > \endverbatim */
637 /* > */
638 /* > \param[in] CAPPLY */
639 /* > \verbatim */
640 /* >          CAPPLY is LOGICAL */
641 /* >     If .TRUE. then access the vector C in the formula above. */
642 /* > \endverbatim */
643 /* > */
644 /* > \param[out] INFO */
645 /* > \verbatim */
646 /* >          INFO is INTEGER */
647 /* >       = 0:  Successful exit. */
648 /* >     i > 0:  The ith argument is invalid. */
649 /* > \endverbatim */
650 /* > */
651 /* > \param[out] WORK */
652 /* > \verbatim */
653 /* >          WORK is COMPLEX array, dimension (2*N). */
654 /* >     Workspace. */
655 /* > \endverbatim */
656 /* > */
657 /* > \param[out] RWORK */
658 /* > \verbatim */
659 /* >          RWORK is REAL array, dimension (N). */
660 /* >     Workspace. */
661 /* > \endverbatim */
662
663 /*  Authors: */
664 /*  ======== */
665
666 /* > \author Univ. of Tennessee */
667 /* > \author Univ. of California Berkeley */
668 /* > \author Univ. of Colorado Denver */
669 /* > \author NAG Ltd. */
670
671 /* > \date December 2016 */
672
673 /* > \ingroup complexGBcomputational */
674
675 /*  ===================================================================== */
676 real cla_gbrcond_c_(char *trans, integer *n, integer *kl, integer *ku, 
677         complex *ab, integer *ldab, complex *afb, integer *ldafb, integer *
678         ipiv, real *c__, logical *capply, integer *info, complex *work, real *
679         rwork)
680 {
681     /* System generated locals */
682     integer ab_dim1, ab_offset, afb_dim1, afb_offset, i__1, i__2, i__3, i__4;
683     real ret_val, r__1, r__2;
684     complex q__1;
685
686     /* Local variables */
687     integer kase, i__, j;
688     extern logical lsame_(char *, char *);
689     integer isave[3];
690     real anorm;
691     extern /* Subroutine */ int clacn2_(integer *, complex *, complex *, real 
692             *, integer *, integer *);
693     integer kd, ke;
694     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen), cgbtrs_(
695             char *, integer *, integer *, integer *, integer *, complex *, 
696             integer *, integer *, complex *, integer *, integer *);
697     real ainvnm, tmp;
698     logical notrans;
699
700
701 /*  -- LAPACK computational routine (version 3.7.0) -- */
702 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
703 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
704 /*     December 2016 */
705
706
707 /*  ===================================================================== */
708
709     /* Parameter adjustments */
710     ab_dim1 = *ldab;
711     ab_offset = 1 + ab_dim1 * 1;
712     ab -= ab_offset;
713     afb_dim1 = *ldafb;
714     afb_offset = 1 + afb_dim1 * 1;
715     afb -= afb_offset;
716     --ipiv;
717     --c__;
718     --work;
719     --rwork;
720
721     /* Function Body */
722     ret_val = 0.f;
723
724     *info = 0;
725     notrans = lsame_(trans, "N");
726     if (! notrans && ! lsame_(trans, "T") && ! lsame_(
727             trans, "C")) {
728         *info = -1;
729     } else if (*n < 0) {
730         *info = -2;
731     } else if (*kl < 0 || *kl > *n - 1) {
732         *info = -3;
733     } else if (*ku < 0 || *ku > *n - 1) {
734         *info = -4;
735     } else if (*ldab < *kl + *ku + 1) {
736         *info = -6;
737     } else if (*ldafb < (*kl << 1) + *ku + 1) {
738         *info = -8;
739     }
740     if (*info != 0) {
741         i__1 = -(*info);
742         xerbla_("CLA_GBRCOND_C", &i__1, (ftnlen)13);
743         return ret_val;
744     }
745
746 /*     Compute norm of op(A)*op2(C). */
747
748     anorm = 0.f;
749     kd = *ku + 1;
750     ke = *kl + 1;
751     if (notrans) {
752         i__1 = *n;
753         for (i__ = 1; i__ <= i__1; ++i__) {
754             tmp = 0.f;
755             if (*capply) {
756 /* Computing MAX */
757                 i__2 = i__ - *kl;
758 /* Computing MIN */
759                 i__4 = i__ + *ku;
760                 i__3 = f2cmin(i__4,*n);
761                 for (j = f2cmax(i__2,1); j <= i__3; ++j) {
762                     i__2 = kd + i__ - j + j * ab_dim1;
763                     tmp += ((r__1 = ab[i__2].r, abs(r__1)) + (r__2 = r_imag(&
764                             ab[kd + i__ - j + j * ab_dim1]), abs(r__2))) / 
765                             c__[j];
766                 }
767             } else {
768 /* Computing MAX */
769                 i__3 = i__ - *kl;
770 /* Computing MIN */
771                 i__4 = i__ + *ku;
772                 i__2 = f2cmin(i__4,*n);
773                 for (j = f2cmax(i__3,1); j <= i__2; ++j) {
774                     i__3 = kd + i__ - j + j * ab_dim1;
775                     tmp += (r__1 = ab[i__3].r, abs(r__1)) + (r__2 = r_imag(&
776                             ab[kd + i__ - j + j * ab_dim1]), abs(r__2));
777                 }
778             }
779             rwork[i__] = tmp;
780             anorm = f2cmax(anorm,tmp);
781         }
782     } else {
783         i__1 = *n;
784         for (i__ = 1; i__ <= i__1; ++i__) {
785             tmp = 0.f;
786             if (*capply) {
787 /* Computing MAX */
788                 i__2 = i__ - *kl;
789 /* Computing MIN */
790                 i__4 = i__ + *ku;
791                 i__3 = f2cmin(i__4,*n);
792                 for (j = f2cmax(i__2,1); j <= i__3; ++j) {
793                     i__2 = ke - i__ + j + i__ * ab_dim1;
794                     tmp += ((r__1 = ab[i__2].r, abs(r__1)) + (r__2 = r_imag(&
795                             ab[ke - i__ + j + i__ * ab_dim1]), abs(r__2))) / 
796                             c__[j];
797                 }
798             } else {
799 /* Computing MAX */
800                 i__3 = i__ - *kl;
801 /* Computing MIN */
802                 i__4 = i__ + *ku;
803                 i__2 = f2cmin(i__4,*n);
804                 for (j = f2cmax(i__3,1); j <= i__2; ++j) {
805                     i__3 = ke - i__ + j + i__ * ab_dim1;
806                     tmp += (r__1 = ab[i__3].r, abs(r__1)) + (r__2 = r_imag(&
807                             ab[ke - i__ + j + i__ * ab_dim1]), abs(r__2));
808                 }
809             }
810             rwork[i__] = tmp;
811             anorm = f2cmax(anorm,tmp);
812         }
813     }
814
815 /*     Quick return if possible. */
816
817     if (*n == 0) {
818         ret_val = 1.f;
819         return ret_val;
820     } else if (anorm == 0.f) {
821         return ret_val;
822     }
823
824 /*     Estimate the norm of inv(op(A)). */
825
826     ainvnm = 0.f;
827
828     kase = 0;
829 L10:
830     clacn2_(n, &work[*n + 1], &work[1], &ainvnm, &kase, isave);
831     if (kase != 0) {
832         if (kase == 2) {
833
834 /*           Multiply by R. */
835
836             i__1 = *n;
837             for (i__ = 1; i__ <= i__1; ++i__) {
838                 i__2 = i__;
839                 i__3 = i__;
840                 i__4 = i__;
841                 q__1.r = rwork[i__4] * work[i__3].r, q__1.i = rwork[i__4] * 
842                         work[i__3].i;
843                 work[i__2].r = q__1.r, work[i__2].i = q__1.i;
844             }
845
846             if (notrans) {
847                 cgbtrs_("No transpose", n, kl, ku, &c__1, &afb[afb_offset], 
848                         ldafb, &ipiv[1], &work[1], n, info);
849             } else {
850                 cgbtrs_("Conjugate transpose", n, kl, ku, &c__1, &afb[
851                         afb_offset], ldafb, &ipiv[1], &work[1], n, info);
852             }
853
854 /*           Multiply by inv(C). */
855
856             if (*capply) {
857                 i__1 = *n;
858                 for (i__ = 1; i__ <= i__1; ++i__) {
859                     i__2 = i__;
860                     i__3 = i__;
861                     i__4 = i__;
862                     q__1.r = c__[i__4] * work[i__3].r, q__1.i = c__[i__4] * 
863                             work[i__3].i;
864                     work[i__2].r = q__1.r, work[i__2].i = q__1.i;
865                 }
866             }
867         } else {
868
869 /*           Multiply by inv(C**H). */
870
871             if (*capply) {
872                 i__1 = *n;
873                 for (i__ = 1; i__ <= i__1; ++i__) {
874                     i__2 = i__;
875                     i__3 = i__;
876                     i__4 = i__;
877                     q__1.r = c__[i__4] * work[i__3].r, q__1.i = c__[i__4] * 
878                             work[i__3].i;
879                     work[i__2].r = q__1.r, work[i__2].i = q__1.i;
880                 }
881             }
882
883             if (notrans) {
884                 cgbtrs_("Conjugate transpose", n, kl, ku, &c__1, &afb[
885                         afb_offset], ldafb, &ipiv[1], &work[1], n, info);
886             } else {
887                 cgbtrs_("No transpose", n, kl, ku, &c__1, &afb[afb_offset], 
888                         ldafb, &ipiv[1], &work[1], n, info);
889             }
890
891 /*           Multiply by R. */
892
893             i__1 = *n;
894             for (i__ = 1; i__ <= i__1; ++i__) {
895                 i__2 = i__;
896                 i__3 = i__;
897                 i__4 = i__;
898                 q__1.r = rwork[i__4] * work[i__3].r, q__1.i = rwork[i__4] * 
899                         work[i__3].i;
900                 work[i__2].r = q__1.r, work[i__2].i = q__1.i;
901             }
902         }
903         goto L10;
904     }
905
906 /*     Compute the estimate of the reciprocal condition number. */
907
908     if (ainvnm != 0.f) {
909         ret_val = 1.f / ainvnm;
910     }
911
912     return ret_val;
913
914 } /* cla_gbrcond_c__ */
915