C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / ctgsna.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 complex c_b19 = {1.f,0.f};
517 static complex c_b20 = {0.f,0.f};
518 static logical c_false = FALSE_;
519 static integer c__3 = 3;
520
521 /* > \brief \b CTGSNA */
522
523 /*  =========== DOCUMENTATION =========== */
524
525 /* Online html documentation available at */
526 /*            http://www.netlib.org/lapack/explore-html/ */
527
528 /* > \htmlonly */
529 /* > Download CTGSNA + dependencies */
530 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ctgsna.
531 f"> */
532 /* > [TGZ]</a> */
533 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ctgsna.
534 f"> */
535 /* > [ZIP]</a> */
536 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ctgsna.
537 f"> */
538 /* > [TXT]</a> */
539 /* > \endhtmlonly */
540
541 /*  Definition: */
542 /*  =========== */
543
544 /*       SUBROUTINE CTGSNA( JOB, HOWMNY, SELECT, N, A, LDA, B, LDB, VL, */
545 /*                          LDVL, VR, LDVR, S, DIF, MM, M, WORK, LWORK, */
546 /*                          IWORK, INFO ) */
547
548 /*       CHARACTER          HOWMNY, JOB */
549 /*       INTEGER            INFO, LDA, LDB, LDVL, LDVR, LWORK, M, MM, N */
550 /*       LOGICAL            SELECT( * ) */
551 /*       INTEGER            IWORK( * ) */
552 /*       REAL               DIF( * ), S( * ) */
553 /*       COMPLEX            A( LDA, * ), B( LDB, * ), VL( LDVL, * ), */
554 /*      $                   VR( LDVR, * ), WORK( * ) */
555
556
557 /* > \par Purpose: */
558 /*  ============= */
559 /* > */
560 /* > \verbatim */
561 /* > */
562 /* > CTGSNA estimates reciprocal condition numbers for specified */
563 /* > eigenvalues and/or eigenvectors of a matrix pair (A, B). */
564 /* > */
565 /* > (A, B) must be in generalized Schur canonical form, that is, A and */
566 /* > B are both upper triangular. */
567 /* > \endverbatim */
568
569 /*  Arguments: */
570 /*  ========== */
571
572 /* > \param[in] JOB */
573 /* > \verbatim */
574 /* >          JOB is CHARACTER*1 */
575 /* >          Specifies whether condition numbers are required for */
576 /* >          eigenvalues (S) or eigenvectors (DIF): */
577 /* >          = 'E': for eigenvalues only (S); */
578 /* >          = 'V': for eigenvectors only (DIF); */
579 /* >          = 'B': for both eigenvalues and eigenvectors (S and DIF). */
580 /* > \endverbatim */
581 /* > */
582 /* > \param[in] HOWMNY */
583 /* > \verbatim */
584 /* >          HOWMNY is CHARACTER*1 */
585 /* >          = 'A': compute condition numbers for all eigenpairs; */
586 /* >          = 'S': compute condition numbers for selected eigenpairs */
587 /* >                 specified by the array SELECT. */
588 /* > \endverbatim */
589 /* > */
590 /* > \param[in] SELECT */
591 /* > \verbatim */
592 /* >          SELECT is LOGICAL array, dimension (N) */
593 /* >          If HOWMNY = 'S', SELECT specifies the eigenpairs for which */
594 /* >          condition numbers are required. To select condition numbers */
595 /* >          for the corresponding j-th eigenvalue and/or eigenvector, */
596 /* >          SELECT(j) must be set to .TRUE.. */
597 /* >          If HOWMNY = 'A', SELECT is not referenced. */
598 /* > \endverbatim */
599 /* > */
600 /* > \param[in] N */
601 /* > \verbatim */
602 /* >          N is INTEGER */
603 /* >          The order of the square matrix pair (A, B). N >= 0. */
604 /* > \endverbatim */
605 /* > */
606 /* > \param[in] A */
607 /* > \verbatim */
608 /* >          A is COMPLEX array, dimension (LDA,N) */
609 /* >          The upper triangular matrix A in the pair (A,B). */
610 /* > \endverbatim */
611 /* > */
612 /* > \param[in] LDA */
613 /* > \verbatim */
614 /* >          LDA is INTEGER */
615 /* >          The leading dimension of the array A. LDA >= f2cmax(1,N). */
616 /* > \endverbatim */
617 /* > */
618 /* > \param[in] B */
619 /* > \verbatim */
620 /* >          B is COMPLEX array, dimension (LDB,N) */
621 /* >          The upper triangular matrix B in the pair (A, B). */
622 /* > \endverbatim */
623 /* > */
624 /* > \param[in] LDB */
625 /* > \verbatim */
626 /* >          LDB is INTEGER */
627 /* >          The leading dimension of the array B. LDB >= f2cmax(1,N). */
628 /* > \endverbatim */
629 /* > */
630 /* > \param[in] VL */
631 /* > \verbatim */
632 /* >          VL is COMPLEX array, dimension (LDVL,M) */
633 /* >          IF JOB = 'E' or 'B', VL must contain left eigenvectors of */
634 /* >          (A, B), corresponding to the eigenpairs specified by HOWMNY */
635 /* >          and SELECT.  The eigenvectors must be stored in consecutive */
636 /* >          columns of VL, as returned by CTGEVC. */
637 /* >          If JOB = 'V', VL is not referenced. */
638 /* > \endverbatim */
639 /* > */
640 /* > \param[in] LDVL */
641 /* > \verbatim */
642 /* >          LDVL is INTEGER */
643 /* >          The leading dimension of the array VL. LDVL >= 1; and */
644 /* >          If JOB = 'E' or 'B', LDVL >= N. */
645 /* > \endverbatim */
646 /* > */
647 /* > \param[in] VR */
648 /* > \verbatim */
649 /* >          VR is COMPLEX array, dimension (LDVR,M) */
650 /* >          IF JOB = 'E' or 'B', VR must contain right eigenvectors of */
651 /* >          (A, B), corresponding to the eigenpairs specified by HOWMNY */
652 /* >          and SELECT.  The eigenvectors must be stored in consecutive */
653 /* >          columns of VR, as returned by CTGEVC. */
654 /* >          If JOB = 'V', VR is not referenced. */
655 /* > \endverbatim */
656 /* > */
657 /* > \param[in] LDVR */
658 /* > \verbatim */
659 /* >          LDVR is INTEGER */
660 /* >          The leading dimension of the array VR. LDVR >= 1; */
661 /* >          If JOB = 'E' or 'B', LDVR >= N. */
662 /* > \endverbatim */
663 /* > */
664 /* > \param[out] S */
665 /* > \verbatim */
666 /* >          S is REAL array, dimension (MM) */
667 /* >          If JOB = 'E' or 'B', the reciprocal condition numbers of the */
668 /* >          selected eigenvalues, stored in consecutive elements of the */
669 /* >          array. */
670 /* >          If JOB = 'V', S is not referenced. */
671 /* > \endverbatim */
672 /* > */
673 /* > \param[out] DIF */
674 /* > \verbatim */
675 /* >          DIF is REAL array, dimension (MM) */
676 /* >          If JOB = 'V' or 'B', the estimated reciprocal condition */
677 /* >          numbers of the selected eigenvectors, stored in consecutive */
678 /* >          elements of the array. */
679 /* >          If the eigenvalues cannot be reordered to compute DIF(j), */
680 /* >          DIF(j) is set to 0; this can only occur when the true value */
681 /* >          would be very small anyway. */
682 /* >          For each eigenvalue/vector specified by SELECT, DIF stores */
683 /* >          a Frobenius norm-based estimate of Difl. */
684 /* >          If JOB = 'E', DIF is not referenced. */
685 /* > \endverbatim */
686 /* > */
687 /* > \param[in] MM */
688 /* > \verbatim */
689 /* >          MM is INTEGER */
690 /* >          The number of elements in the arrays S and DIF. MM >= M. */
691 /* > \endverbatim */
692 /* > */
693 /* > \param[out] M */
694 /* > \verbatim */
695 /* >          M is INTEGER */
696 /* >          The number of elements of the arrays S and DIF used to store */
697 /* >          the specified condition numbers; for each selected eigenvalue */
698 /* >          one element is used. If HOWMNY = 'A', M is set to N. */
699 /* > \endverbatim */
700 /* > */
701 /* > \param[out] WORK */
702 /* > \verbatim */
703 /* >          WORK is COMPLEX array, dimension (MAX(1,LWORK)) */
704 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
705 /* > \endverbatim */
706 /* > */
707 /* > \param[in] LWORK */
708 /* > \verbatim */
709 /* >          LWORK is INTEGER */
710 /* >          The dimension of the array WORK. LWORK >= f2cmax(1,N). */
711 /* >          If JOB = 'V' or 'B', LWORK >= f2cmax(1,2*N*N). */
712 /* > \endverbatim */
713 /* > */
714 /* > \param[out] IWORK */
715 /* > \verbatim */
716 /* >          IWORK is INTEGER array, dimension (N+2) */
717 /* >          If JOB = 'E', IWORK is not referenced. */
718 /* > \endverbatim */
719 /* > */
720 /* > \param[out] INFO */
721 /* > \verbatim */
722 /* >          INFO is INTEGER */
723 /* >          = 0: Successful exit */
724 /* >          < 0: If INFO = -i, the i-th argument had an illegal value */
725 /* > \endverbatim */
726
727 /*  Authors: */
728 /*  ======== */
729
730 /* > \author Univ. of Tennessee */
731 /* > \author Univ. of California Berkeley */
732 /* > \author Univ. of Colorado Denver */
733 /* > \author NAG Ltd. */
734
735 /* > \date December 2016 */
736
737 /* > \ingroup complexOTHERcomputational */
738
739 /* > \par Further Details: */
740 /*  ===================== */
741 /* > */
742 /* > \verbatim */
743 /* > */
744 /* >  The reciprocal of the condition number of the i-th generalized */
745 /* >  eigenvalue w = (a, b) is defined as */
746 /* > */
747 /* >          S(I) = (|v**HAu|**2 + |v**HBu|**2)**(1/2) / (norm(u)*norm(v)) */
748 /* > */
749 /* >  where u and v are the right and left eigenvectors of (A, B) */
750 /* >  corresponding to w; |z| denotes the absolute value of the complex */
751 /* >  number, and norm(u) denotes the 2-norm of the vector u. The pair */
752 /* >  (a, b) corresponds to an eigenvalue w = a/b (= v**HAu/v**HBu) of the */
753 /* >  matrix pair (A, B). If both a and b equal zero, then (A,B) is */
754 /* >  singular and S(I) = -1 is returned. */
755 /* > */
756 /* >  An approximate error bound on the chordal distance between the i-th */
757 /* >  computed generalized eigenvalue w and the corresponding exact */
758 /* >  eigenvalue lambda is */
759 /* > */
760 /* >          chord(w, lambda) <=   EPS * norm(A, B) / S(I), */
761 /* > */
762 /* >  where EPS is the machine precision. */
763 /* > */
764 /* >  The reciprocal of the condition number of the right eigenvector u */
765 /* >  and left eigenvector v corresponding to the generalized eigenvalue w */
766 /* >  is defined as follows. Suppose */
767 /* > */
768 /* >                   (A, B) = ( a   *  ) ( b  *  )  1 */
769 /* >                            ( 0  A22 ),( 0 B22 )  n-1 */
770 /* >                              1  n-1     1 n-1 */
771 /* > */
772 /* >  Then the reciprocal condition number DIF(I) is */
773 /* > */
774 /* >          Difl[(a, b), (A22, B22)]  = sigma-f2cmin( Zl ) */
775 /* > */
776 /* >  where sigma-f2cmin(Zl) denotes the smallest singular value of */
777 /* > */
778 /* >         Zl = [ kron(a, In-1) -kron(1, A22) ] */
779 /* >              [ kron(b, In-1) -kron(1, B22) ]. */
780 /* > */
781 /* >  Here In-1 is the identity matrix of size n-1 and X**H is the conjugate */
782 /* >  transpose of X. kron(X, Y) is the Kronecker product between the */
783 /* >  matrices X and Y. */
784 /* > */
785 /* >  We approximate the smallest singular value of Zl with an upper */
786 /* >  bound. This is done by CLATDF. */
787 /* > */
788 /* >  An approximate error bound for a computed eigenvector VL(i) or */
789 /* >  VR(i) is given by */
790 /* > */
791 /* >                      EPS * norm(A, B) / DIF(i). */
792 /* > */
793 /* >  See ref. [2-3] for more details and further references. */
794 /* > \endverbatim */
795
796 /* > \par Contributors: */
797 /*  ================== */
798 /* > */
799 /* >     Bo Kagstrom and Peter Poromaa, Department of Computing Science, */
800 /* >     Umea University, S-901 87 Umea, Sweden. */
801
802 /* > \par References: */
803 /*  ================ */
804 /* > */
805 /* > \verbatim */
806 /* > */
807 /* >  [1] B. Kagstrom; A Direct Method for Reordering Eigenvalues in the */
808 /* >      Generalized Real Schur Form of a Regular Matrix Pair (A, B), in */
809 /* >      M.S. Moonen et al (eds), Linear Algebra for Large Scale and */
810 /* >      Real-Time Applications, Kluwer Academic Publ. 1993, pp 195-218. */
811 /* > */
812 /* >  [2] B. Kagstrom and P. Poromaa; Computing Eigenspaces with Specified */
813 /* >      Eigenvalues of a Regular Matrix Pair (A, B) and Condition */
814 /* >      Estimation: Theory, Algorithms and Software, Report */
815 /* >      UMINF - 94.04, Department of Computing Science, Umea University, */
816 /* >      S-901 87 Umea, Sweden, 1994. Also as LAPACK Working Note 87. */
817 /* >      To appear in Numerical Algorithms, 1996. */
818 /* > */
819 /* >  [3] B. Kagstrom and P. Poromaa, LAPACK-Style Algorithms and Software */
820 /* >      for Solving the Generalized Sylvester Equation and Estimating the */
821 /* >      Separation between Regular Matrix Pairs, Report UMINF - 93.23, */
822 /* >      Department of Computing Science, Umea University, S-901 87 Umea, */
823 /* >      Sweden, December 1993, Revised April 1994, Also as LAPACK Working */
824 /* >      Note 75. */
825 /* >      To appear in ACM Trans. on Math. Software, Vol 22, No 1, 1996. */
826 /* > \endverbatim */
827 /* > */
828 /*  ===================================================================== */
829 /* Subroutine */ int ctgsna_(char *job, char *howmny, logical *select, 
830         integer *n, complex *a, integer *lda, complex *b, integer *ldb, 
831         complex *vl, integer *ldvl, complex *vr, integer *ldvr, real *s, real 
832         *dif, integer *mm, integer *m, complex *work, integer *lwork, integer 
833         *iwork, integer *info)
834 {
835     /* System generated locals */
836     integer a_dim1, a_offset, b_dim1, b_offset, vl_dim1, vl_offset, vr_dim1, 
837             vr_offset, i__1;
838     real r__1, r__2;
839     complex q__1;
840
841     /* Local variables */
842     real cond;
843     integer ierr, ifst;
844     real lnrm;
845     complex yhax, yhbx;
846     integer ilst;
847     real rnrm;
848     integer i__, k;
849     real scale;
850     extern /* Complex */ VOID cdotc_(complex *, integer *, complex *, integer 
851             *, complex *, integer *);
852     extern logical lsame_(char *, char *);
853     extern /* Subroutine */ int cgemv_(char *, integer *, integer *, complex *
854             , complex *, integer *, complex *, integer *, complex *, complex *
855             , integer *);
856     integer lwmin;
857     logical wants;
858     complex dummy[1];
859     integer n1, n2;
860     extern real scnrm2_(integer *, complex *, integer *), slapy2_(real *, 
861             real *);
862     complex dummy1[1];
863     extern /* Subroutine */ int slabad_(real *, real *);
864     integer ks;
865     extern real slamch_(char *);
866     extern /* Subroutine */ int clacpy_(char *, integer *, integer *, complex 
867             *, integer *, complex *, integer *), ctgexc_(logical *, 
868             logical *, integer *, complex *, integer *, complex *, integer *, 
869             complex *, integer *, complex *, integer *, integer *, integer *, 
870             integer *), xerbla_(char *, integer *, ftnlen);
871     real bignum;
872     logical wantbh, wantdf, somcon;
873     extern /* Subroutine */ int ctgsyl_(char *, integer *, integer *, integer 
874             *, complex *, integer *, complex *, integer *, complex *, integer 
875             *, complex *, integer *, complex *, integer *, complex *, integer 
876             *, real *, real *, complex *, integer *, integer *, integer *);
877     real smlnum;
878     logical lquery;
879     real eps;
880
881
882 /*  -- LAPACK computational routine (version 3.7.0) -- */
883 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
884 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
885 /*     December 2016 */
886
887
888 /*  ===================================================================== */
889
890
891 /*     Decode and test the input parameters */
892
893     /* Parameter adjustments */
894     --select;
895     a_dim1 = *lda;
896     a_offset = 1 + a_dim1 * 1;
897     a -= a_offset;
898     b_dim1 = *ldb;
899     b_offset = 1 + b_dim1 * 1;
900     b -= b_offset;
901     vl_dim1 = *ldvl;
902     vl_offset = 1 + vl_dim1 * 1;
903     vl -= vl_offset;
904     vr_dim1 = *ldvr;
905     vr_offset = 1 + vr_dim1 * 1;
906     vr -= vr_offset;
907     --s;
908     --dif;
909     --work;
910     --iwork;
911
912     /* Function Body */
913     wantbh = lsame_(job, "B");
914     wants = lsame_(job, "E") || wantbh;
915     wantdf = lsame_(job, "V") || wantbh;
916
917     somcon = lsame_(howmny, "S");
918
919     *info = 0;
920     lquery = *lwork == -1;
921
922     if (! wants && ! wantdf) {
923         *info = -1;
924     } else if (! lsame_(howmny, "A") && ! somcon) {
925         *info = -2;
926     } else if (*n < 0) {
927         *info = -4;
928     } else if (*lda < f2cmax(1,*n)) {
929         *info = -6;
930     } else if (*ldb < f2cmax(1,*n)) {
931         *info = -8;
932     } else if (wants && *ldvl < *n) {
933         *info = -10;
934     } else if (wants && *ldvr < *n) {
935         *info = -12;
936     } else {
937
938 /*        Set M to the number of eigenpairs for which condition numbers */
939 /*        are required, and test MM. */
940
941         if (somcon) {
942             *m = 0;
943             i__1 = *n;
944             for (k = 1; k <= i__1; ++k) {
945                 if (select[k]) {
946                     ++(*m);
947                 }
948 /* L10: */
949             }
950         } else {
951             *m = *n;
952         }
953
954         if (*n == 0) {
955             lwmin = 1;
956         } else if (lsame_(job, "V") || lsame_(job, 
957                 "B")) {
958             lwmin = (*n << 1) * *n;
959         } else {
960             lwmin = *n;
961         }
962         work[1].r = (real) lwmin, work[1].i = 0.f;
963
964         if (*mm < *m) {
965             *info = -15;
966         } else if (*lwork < lwmin && ! lquery) {
967             *info = -18;
968         }
969     }
970
971     if (*info != 0) {
972         i__1 = -(*info);
973         xerbla_("CTGSNA", &i__1, (ftnlen)6);
974         return 0;
975     } else if (lquery) {
976         return 0;
977     }
978
979 /*     Quick return if possible */
980
981     if (*n == 0) {
982         return 0;
983     }
984
985 /*     Get machine constants */
986
987     eps = slamch_("P");
988     smlnum = slamch_("S") / eps;
989     bignum = 1.f / smlnum;
990     slabad_(&smlnum, &bignum);
991     ks = 0;
992     i__1 = *n;
993     for (k = 1; k <= i__1; ++k) {
994
995 /*        Determine whether condition numbers are required for the k-th */
996 /*        eigenpair. */
997
998         if (somcon) {
999             if (! select[k]) {
1000                 goto L20;
1001             }
1002         }
1003
1004         ++ks;
1005
1006         if (wants) {
1007
1008 /*           Compute the reciprocal condition number of the k-th */
1009 /*           eigenvalue. */
1010
1011             rnrm = scnrm2_(n, &vr[ks * vr_dim1 + 1], &c__1);
1012             lnrm = scnrm2_(n, &vl[ks * vl_dim1 + 1], &c__1);
1013             cgemv_("N", n, n, &c_b19, &a[a_offset], lda, &vr[ks * vr_dim1 + 1]
1014                     , &c__1, &c_b20, &work[1], &c__1);
1015             cdotc_(&q__1, n, &work[1], &c__1, &vl[ks * vl_dim1 + 1], &c__1);
1016             yhax.r = q__1.r, yhax.i = q__1.i;
1017             cgemv_("N", n, n, &c_b19, &b[b_offset], ldb, &vr[ks * vr_dim1 + 1]
1018                     , &c__1, &c_b20, &work[1], &c__1);
1019             cdotc_(&q__1, n, &work[1], &c__1, &vl[ks * vl_dim1 + 1], &c__1);
1020             yhbx.r = q__1.r, yhbx.i = q__1.i;
1021             r__1 = c_abs(&yhax);
1022             r__2 = c_abs(&yhbx);
1023             cond = slapy2_(&r__1, &r__2);
1024             if (cond == 0.f) {
1025                 s[ks] = -1.f;
1026             } else {
1027                 s[ks] = cond / (rnrm * lnrm);
1028             }
1029         }
1030
1031         if (wantdf) {
1032             if (*n == 1) {
1033                 r__1 = c_abs(&a[a_dim1 + 1]);
1034                 r__2 = c_abs(&b[b_dim1 + 1]);
1035                 dif[ks] = slapy2_(&r__1, &r__2);
1036             } else {
1037
1038 /*              Estimate the reciprocal condition number of the k-th */
1039 /*              eigenvectors. */
1040
1041 /*              Copy the matrix (A, B) to the array WORK and move the */
1042 /*              (k,k)th pair to the (1,1) position. */
1043
1044                 clacpy_("Full", n, n, &a[a_offset], lda, &work[1], n);
1045                 clacpy_("Full", n, n, &b[b_offset], ldb, &work[*n * *n + 1], 
1046                         n);
1047                 ifst = k;
1048                 ilst = 1;
1049
1050                 ctgexc_(&c_false, &c_false, n, &work[1], n, &work[*n * *n + 1]
1051                         , n, dummy, &c__1, dummy1, &c__1, &ifst, &ilst, &ierr)
1052                         ;
1053
1054                 if (ierr > 0) {
1055
1056 /*                 Ill-conditioned problem - swap rejected. */
1057
1058                     dif[ks] = 0.f;
1059                 } else {
1060
1061 /*                 Reordering successful, solve generalized Sylvester */
1062 /*                 equation for R and L, */
1063 /*                            A22 * R - L * A11 = A12 */
1064 /*                            B22 * R - L * B11 = B12, */
1065 /*                 and compute estimate of Difl[(A11,B11), (A22, B22)]. */
1066
1067                     n1 = 1;
1068                     n2 = *n - n1;
1069                     i__ = *n * *n + 1;
1070                     ctgsyl_("N", &c__3, &n2, &n1, &work[*n * n1 + n1 + 1], n, 
1071                             &work[1], n, &work[n1 + 1], n, &work[*n * n1 + n1 
1072                             + i__], n, &work[i__], n, &work[n1 + i__], n, &
1073                             scale, &dif[ks], dummy, &c__1, &iwork[1], &ierr);
1074                 }
1075             }
1076         }
1077
1078 L20:
1079         ;
1080     }
1081     work[1].r = (real) lwmin, work[1].i = 0.f;
1082     return 0;
1083
1084 /*     End of CTGSNA */
1085
1086 } /* ctgsna_ */
1087