C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / cggevx.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 complex c_b1 = {0.f,0.f};
516 static complex c_b2 = {1.f,0.f};
517 static integer c__1 = 1;
518 static integer c__0 = 0;
519
520 /* > \brief <b> CGGEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE mat
521 rices</b> */
522
523 /*  =========== DOCUMENTATION =========== */
524
525 /* Online html documentation available at */
526 /*            http://www.netlib.org/lapack/explore-html/ */
527
528 /* > \htmlonly */
529 /* > Download CGGEVX + dependencies */
530 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cggevx.
531 f"> */
532 /* > [TGZ]</a> */
533 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cggevx.
534 f"> */
535 /* > [ZIP]</a> */
536 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cggevx.
537 f"> */
538 /* > [TXT]</a> */
539 /* > \endhtmlonly */
540
541 /*  Definition: */
542 /*  =========== */
543
544 /*       SUBROUTINE CGGEVX( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB, */
545 /*                          ALPHA, BETA, VL, LDVL, VR, LDVR, ILO, IHI, */
546 /*                          LSCALE, RSCALE, ABNRM, BBNRM, RCONDE, RCONDV, */
547 /*                          WORK, LWORK, RWORK, IWORK, BWORK, INFO ) */
548
549 /*       CHARACTER          BALANC, JOBVL, JOBVR, SENSE */
550 /*       INTEGER            IHI, ILO, INFO, LDA, LDB, LDVL, LDVR, LWORK, N */
551 /*       REAL               ABNRM, BBNRM */
552 /*       LOGICAL            BWORK( * ) */
553 /*       INTEGER            IWORK( * ) */
554 /*       REAL               LSCALE( * ), RCONDE( * ), RCONDV( * ), */
555 /*      $                   RSCALE( * ), RWORK( * ) */
556 /*       COMPLEX            A( LDA, * ), ALPHA( * ), B( LDB, * ), */
557 /*      $                   BETA( * ), VL( LDVL, * ), VR( LDVR, * ), */
558 /*      $                   WORK( * ) */
559
560
561 /* > \par Purpose: */
562 /*  ============= */
563 /* > */
564 /* > \verbatim */
565 /* > */
566 /* > CGGEVX computes for a pair of N-by-N complex nonsymmetric matrices */
567 /* > (A,B) the generalized eigenvalues, and optionally, the left and/or */
568 /* > right generalized eigenvectors. */
569 /* > */
570 /* > Optionally, it also computes a balancing transformation to improve */
571 /* > the conditioning of the eigenvalues and eigenvectors (ILO, IHI, */
572 /* > LSCALE, RSCALE, ABNRM, and BBNRM), reciprocal condition numbers for */
573 /* > the eigenvalues (RCONDE), and reciprocal condition numbers for the */
574 /* > right eigenvectors (RCONDV). */
575 /* > */
576 /* > A generalized eigenvalue for a pair of matrices (A,B) is a scalar */
577 /* > lambda or a ratio alpha/beta = lambda, such that A - lambda*B is */
578 /* > singular. It is usually represented as the pair (alpha,beta), as */
579 /* > there is a reasonable interpretation for beta=0, and even for both */
580 /* > being zero. */
581 /* > */
582 /* > The right eigenvector v(j) corresponding to the eigenvalue lambda(j) */
583 /* > of (A,B) satisfies */
584 /* >                  A * v(j) = lambda(j) * B * v(j) . */
585 /* > The left eigenvector u(j) corresponding to the eigenvalue lambda(j) */
586 /* > of (A,B) satisfies */
587 /* >                  u(j)**H * A  = lambda(j) * u(j)**H * B. */
588 /* > where u(j)**H is the conjugate-transpose of u(j). */
589 /* > */
590 /* > \endverbatim */
591
592 /*  Arguments: */
593 /*  ========== */
594
595 /* > \param[in] BALANC */
596 /* > \verbatim */
597 /* >          BALANC is CHARACTER*1 */
598 /* >          Specifies the balance option to be performed: */
599 /* >          = 'N':  do not diagonally scale or permute; */
600 /* >          = 'P':  permute only; */
601 /* >          = 'S':  scale only; */
602 /* >          = 'B':  both permute and scale. */
603 /* >          Computed reciprocal condition numbers will be for the */
604 /* >          matrices after permuting and/or balancing. Permuting does */
605 /* >          not change condition numbers (in exact arithmetic), but */
606 /* >          balancing does. */
607 /* > \endverbatim */
608 /* > */
609 /* > \param[in] JOBVL */
610 /* > \verbatim */
611 /* >          JOBVL is CHARACTER*1 */
612 /* >          = 'N':  do not compute the left generalized eigenvectors; */
613 /* >          = 'V':  compute the left generalized eigenvectors. */
614 /* > \endverbatim */
615 /* > */
616 /* > \param[in] JOBVR */
617 /* > \verbatim */
618 /* >          JOBVR is CHARACTER*1 */
619 /* >          = 'N':  do not compute the right generalized eigenvectors; */
620 /* >          = 'V':  compute the right generalized eigenvectors. */
621 /* > \endverbatim */
622 /* > */
623 /* > \param[in] SENSE */
624 /* > \verbatim */
625 /* >          SENSE is CHARACTER*1 */
626 /* >          Determines which reciprocal condition numbers are computed. */
627 /* >          = 'N': none are computed; */
628 /* >          = 'E': computed for eigenvalues only; */
629 /* >          = 'V': computed for eigenvectors only; */
630 /* >          = 'B': computed for eigenvalues and eigenvectors. */
631 /* > \endverbatim */
632 /* > */
633 /* > \param[in] N */
634 /* > \verbatim */
635 /* >          N is INTEGER */
636 /* >          The order of the matrices A, B, VL, and VR.  N >= 0. */
637 /* > \endverbatim */
638 /* > */
639 /* > \param[in,out] A */
640 /* > \verbatim */
641 /* >          A is COMPLEX array, dimension (LDA, N) */
642 /* >          On entry, the matrix A in the pair (A,B). */
643 /* >          On exit, A has been overwritten. If JOBVL='V' or JOBVR='V' */
644 /* >          or both, then A contains the first part of the complex Schur */
645 /* >          form of the "balanced" versions of the input A and B. */
646 /* > \endverbatim */
647 /* > */
648 /* > \param[in] LDA */
649 /* > \verbatim */
650 /* >          LDA is INTEGER */
651 /* >          The leading dimension of A.  LDA >= f2cmax(1,N). */
652 /* > \endverbatim */
653 /* > */
654 /* > \param[in,out] B */
655 /* > \verbatim */
656 /* >          B is COMPLEX array, dimension (LDB, N) */
657 /* >          On entry, the matrix B in the pair (A,B). */
658 /* >          On exit, B has been overwritten. If JOBVL='V' or JOBVR='V' */
659 /* >          or both, then B contains the second part of the complex */
660 /* >          Schur form of the "balanced" versions of the input A and B. */
661 /* > \endverbatim */
662 /* > */
663 /* > \param[in] LDB */
664 /* > \verbatim */
665 /* >          LDB is INTEGER */
666 /* >          The leading dimension of B.  LDB >= f2cmax(1,N). */
667 /* > \endverbatim */
668 /* > */
669 /* > \param[out] ALPHA */
670 /* > \verbatim */
671 /* >          ALPHA is COMPLEX array, dimension (N) */
672 /* > \endverbatim */
673 /* > */
674 /* > \param[out] BETA */
675 /* > \verbatim */
676 /* >          BETA is COMPLEX array, dimension (N) */
677 /* >          On exit, ALPHA(j)/BETA(j), j=1,...,N, will be the generalized */
678 /* >          eigenvalues. */
679 /* > */
680 /* >          Note: the quotient ALPHA(j)/BETA(j) ) may easily over- or */
681 /* >          underflow, and BETA(j) may even be zero.  Thus, the user */
682 /* >          should avoid naively computing the ratio ALPHA/BETA. */
683 /* >          However, ALPHA will be always less than and usually */
684 /* >          comparable with norm(A) in magnitude, and BETA always less */
685 /* >          than and usually comparable with norm(B). */
686 /* > \endverbatim */
687 /* > */
688 /* > \param[out] VL */
689 /* > \verbatim */
690 /* >          VL is COMPLEX array, dimension (LDVL,N) */
691 /* >          If JOBVL = 'V', the left generalized eigenvectors u(j) are */
692 /* >          stored one after another in the columns of VL, in the same */
693 /* >          order as their eigenvalues. */
694 /* >          Each eigenvector will be scaled so the largest component */
695 /* >          will have abs(real part) + abs(imag. part) = 1. */
696 /* >          Not referenced if JOBVL = 'N'. */
697 /* > \endverbatim */
698 /* > */
699 /* > \param[in] LDVL */
700 /* > \verbatim */
701 /* >          LDVL is INTEGER */
702 /* >          The leading dimension of the matrix VL. LDVL >= 1, and */
703 /* >          if JOBVL = 'V', LDVL >= N. */
704 /* > \endverbatim */
705 /* > */
706 /* > \param[out] VR */
707 /* > \verbatim */
708 /* >          VR is COMPLEX array, dimension (LDVR,N) */
709 /* >          If JOBVR = 'V', the right generalized eigenvectors v(j) are */
710 /* >          stored one after another in the columns of VR, in the same */
711 /* >          order as their eigenvalues. */
712 /* >          Each eigenvector will be scaled so the largest component */
713 /* >          will have abs(real part) + abs(imag. part) = 1. */
714 /* >          Not referenced if JOBVR = 'N'. */
715 /* > \endverbatim */
716 /* > */
717 /* > \param[in] LDVR */
718 /* > \verbatim */
719 /* >          LDVR is INTEGER */
720 /* >          The leading dimension of the matrix VR. LDVR >= 1, and */
721 /* >          if JOBVR = 'V', LDVR >= N. */
722 /* > \endverbatim */
723 /* > */
724 /* > \param[out] ILO */
725 /* > \verbatim */
726 /* >          ILO is INTEGER */
727 /* > \endverbatim */
728 /* > */
729 /* > \param[out] IHI */
730 /* > \verbatim */
731 /* >          IHI is INTEGER */
732 /* >          ILO and IHI are integer values such that on exit */
733 /* >          A(i,j) = 0 and B(i,j) = 0 if i > j and */
734 /* >          j = 1,...,ILO-1 or i = IHI+1,...,N. */
735 /* >          If BALANC = 'N' or 'S', ILO = 1 and IHI = N. */
736 /* > \endverbatim */
737 /* > */
738 /* > \param[out] LSCALE */
739 /* > \verbatim */
740 /* >          LSCALE is REAL array, dimension (N) */
741 /* >          Details of the permutations and scaling factors applied */
742 /* >          to the left side of A and B.  If PL(j) is the index of the */
743 /* >          row interchanged with row j, and DL(j) is the scaling */
744 /* >          factor applied to row j, then */
745 /* >            LSCALE(j) = PL(j)  for j = 1,...,ILO-1 */
746 /* >                      = DL(j)  for j = ILO,...,IHI */
747 /* >                      = PL(j)  for j = IHI+1,...,N. */
748 /* >          The order in which the interchanges are made is N to IHI+1, */
749 /* >          then 1 to ILO-1. */
750 /* > \endverbatim */
751 /* > */
752 /* > \param[out] RSCALE */
753 /* > \verbatim */
754 /* >          RSCALE is REAL array, dimension (N) */
755 /* >          Details of the permutations and scaling factors applied */
756 /* >          to the right side of A and B.  If PR(j) is the index of the */
757 /* >          column interchanged with column j, and DR(j) is the scaling */
758 /* >          factor applied to column j, then */
759 /* >            RSCALE(j) = PR(j)  for j = 1,...,ILO-1 */
760 /* >                      = DR(j)  for j = ILO,...,IHI */
761 /* >                      = PR(j)  for j = IHI+1,...,N */
762 /* >          The order in which the interchanges are made is N to IHI+1, */
763 /* >          then 1 to ILO-1. */
764 /* > \endverbatim */
765 /* > */
766 /* > \param[out] ABNRM */
767 /* > \verbatim */
768 /* >          ABNRM is REAL */
769 /* >          The one-norm of the balanced matrix A. */
770 /* > \endverbatim */
771 /* > */
772 /* > \param[out] BBNRM */
773 /* > \verbatim */
774 /* >          BBNRM is REAL */
775 /* >          The one-norm of the balanced matrix B. */
776 /* > \endverbatim */
777 /* > */
778 /* > \param[out] RCONDE */
779 /* > \verbatim */
780 /* >          RCONDE is REAL array, dimension (N) */
781 /* >          If SENSE = 'E' or 'B', the reciprocal condition numbers of */
782 /* >          the eigenvalues, stored in consecutive elements of the array. */
783 /* >          If SENSE = 'N' or 'V', RCONDE is not referenced. */
784 /* > \endverbatim */
785 /* > */
786 /* > \param[out] RCONDV */
787 /* > \verbatim */
788 /* >          RCONDV is REAL array, dimension (N) */
789 /* >          If SENSE = 'V' or 'B', the estimated reciprocal condition */
790 /* >          numbers of the eigenvectors, stored in consecutive elements */
791 /* >          of the array. If the eigenvalues cannot be reordered to */
792 /* >          compute RCONDV(j), RCONDV(j) is set to 0; this can only occur */
793 /* >          when the true value would be very small anyway. */
794 /* >          If SENSE = 'N' or 'E', RCONDV is not referenced. */
795 /* > \endverbatim */
796 /* > */
797 /* > \param[out] WORK */
798 /* > \verbatim */
799 /* >          WORK is COMPLEX array, dimension (MAX(1,LWORK)) */
800 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
801 /* > \endverbatim */
802 /* > */
803 /* > \param[in] LWORK */
804 /* > \verbatim */
805 /* >          LWORK is INTEGER */
806 /* >          The dimension of the array WORK. LWORK >= f2cmax(1,2*N). */
807 /* >          If SENSE = 'E', LWORK >= f2cmax(1,4*N). */
808 /* >          If SENSE = 'V' or 'B', LWORK >= f2cmax(1,2*N*N+2*N). */
809 /* > */
810 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
811 /* >          only calculates the optimal size of the WORK array, returns */
812 /* >          this value as the first entry of the WORK array, and no error */
813 /* >          message related to LWORK is issued by XERBLA. */
814 /* > \endverbatim */
815 /* > */
816 /* > \param[out] RWORK */
817 /* > \verbatim */
818 /* >          RWORK is REAL array, dimension (lrwork) */
819 /* >          lrwork must be at least f2cmax(1,6*N) if BALANC = 'S' or 'B', */
820 /* >          and at least f2cmax(1,2*N) otherwise. */
821 /* >          Real workspace. */
822 /* > \endverbatim */
823 /* > */
824 /* > \param[out] IWORK */
825 /* > \verbatim */
826 /* >          IWORK is INTEGER array, dimension (N+2) */
827 /* >          If SENSE = 'E', IWORK is not referenced. */
828 /* > \endverbatim */
829 /* > */
830 /* > \param[out] BWORK */
831 /* > \verbatim */
832 /* >          BWORK is LOGICAL array, dimension (N) */
833 /* >          If SENSE = 'N', BWORK is not referenced. */
834 /* > \endverbatim */
835 /* > */
836 /* > \param[out] INFO */
837 /* > \verbatim */
838 /* >          INFO is INTEGER */
839 /* >          = 0:  successful exit */
840 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value. */
841 /* >          = 1,...,N: */
842 /* >                The QZ iteration failed.  No eigenvectors have been */
843 /* >                calculated, but ALPHA(j) and BETA(j) should be correct */
844 /* >                for j=INFO+1,...,N. */
845 /* >          > N:  =N+1: other than QZ iteration failed in CHGEQZ. */
846 /* >                =N+2: error return from CTGEVC. */
847 /* > \endverbatim */
848
849 /*  Authors: */
850 /*  ======== */
851
852 /* > \author Univ. of Tennessee */
853 /* > \author Univ. of California Berkeley */
854 /* > \author Univ. of Colorado Denver */
855 /* > \author NAG Ltd. */
856
857 /* > \date April 2012 */
858
859 /* > \ingroup complexGEeigen */
860
861 /* > \par Further Details: */
862 /*  ===================== */
863 /* > */
864 /* > \verbatim */
865 /* > */
866 /* >  Balancing a matrix pair (A,B) includes, first, permuting rows and */
867 /* >  columns to isolate eigenvalues, second, applying diagonal similarity */
868 /* >  transformation to the rows and columns to make the rows and columns */
869 /* >  as close in norm as possible. The computed reciprocal condition */
870 /* >  numbers correspond to the balanced matrix. Permuting rows and columns */
871 /* >  will not change the condition numbers (in exact arithmetic) but */
872 /* >  diagonal scaling will.  For further explanation of balancing, see */
873 /* >  section 4.11.1.2 of LAPACK Users' Guide. */
874 /* > */
875 /* >  An approximate error bound on the chordal distance between the i-th */
876 /* >  computed generalized eigenvalue w and the corresponding exact */
877 /* >  eigenvalue lambda is */
878 /* > */
879 /* >       chord(w, lambda) <= EPS * norm(ABNRM, BBNRM) / RCONDE(I) */
880 /* > */
881 /* >  An approximate error bound for the angle between the i-th computed */
882 /* >  eigenvector VL(i) or VR(i) is given by */
883 /* > */
884 /* >       EPS * norm(ABNRM, BBNRM) / DIF(i). */
885 /* > */
886 /* >  For further explanation of the reciprocal condition numbers RCONDE */
887 /* >  and RCONDV, see section 4.11 of LAPACK User's Guide. */
888 /* > \endverbatim */
889 /* > */
890 /*  ===================================================================== */
891 /* Subroutine */ int cggevx_(char *balanc, char *jobvl, char *jobvr, char *
892         sense, integer *n, complex *a, integer *lda, complex *b, integer *ldb,
893          complex *alpha, complex *beta, complex *vl, integer *ldvl, complex *
894         vr, integer *ldvr, integer *ilo, integer *ihi, real *lscale, real *
895         rscale, real *abnrm, real *bbnrm, real *rconde, real *rcondv, complex 
896         *work, integer *lwork, real *rwork, integer *iwork, logical *bwork, 
897         integer *info)
898 {
899     /* System generated locals */
900     integer a_dim1, a_offset, b_dim1, b_offset, vl_dim1, vl_offset, vr_dim1, 
901             vr_offset, i__1, i__2, i__3, i__4;
902     real r__1, r__2, r__3, r__4;
903     complex q__1;
904
905     /* Local variables */
906     real anrm, bnrm;
907     integer ierr, itau;
908     real temp;
909     logical ilvl, ilvr;
910     integer iwrk, iwrk1, i__, j, m;
911     extern logical lsame_(char *, char *);
912     integer icols;
913     logical noscl;
914     integer irows, jc;
915     extern /* Subroutine */ int cggbak_(char *, char *, integer *, integer *, 
916             integer *, real *, real *, integer *, complex *, integer *, 
917             integer *), cggbal_(char *, integer *, complex *, 
918             integer *, complex *, integer *, integer *, integer *, real *, 
919             real *, real *, integer *), slabad_(real *, real *);
920     integer in;
921     extern real clange_(char *, integer *, integer *, complex *, integer *, 
922             real *);
923     integer jr;
924     extern /* Subroutine */ int cgghrd_(char *, char *, integer *, integer *, 
925             integer *, complex *, integer *, complex *, integer *, complex *, 
926             integer *, complex *, integer *, integer *), 
927             clascl_(char *, integer *, integer *, real *, real *, integer *, 
928             integer *, complex *, integer *, integer *);
929     logical ilascl, ilbscl;
930     extern /* Subroutine */ int cgeqrf_(integer *, integer *, complex *, 
931             integer *, complex *, complex *, integer *, integer *), clacpy_(
932             char *, integer *, integer *, complex *, integer *, complex *, 
933             integer *), claset_(char *, integer *, integer *, complex 
934             *, complex *, complex *, integer *);
935     logical ldumma[1];
936     char chtemp[1];
937     real bignum;
938     extern /* Subroutine */ int chgeqz_(char *, char *, char *, integer *, 
939             integer *, integer *, complex *, integer *, complex *, integer *, 
940             complex *, complex *, complex *, integer *, complex *, integer *, 
941             complex *, integer *, real *, integer *), 
942             ctgevc_(char *, char *, logical *, integer *, complex *, integer *
943             , complex *, integer *, complex *, integer *, complex *, integer *
944             , integer *, integer *, complex *, real *, integer *);
945     integer ijobvl;
946     extern /* Subroutine */ int ctgsna_(char *, char *, logical *, integer *, 
947             complex *, integer *, complex *, integer *, complex *, integer *, 
948             complex *, integer *, real *, real *, integer *, integer *, 
949             complex *, integer *, integer *, integer *), 
950             slascl_(char *, integer *, integer *, real *, real *, integer *, 
951             integer *, real *, integer *, integer *), xerbla_(char *, 
952             integer *, ftnlen);
953     extern integer ilaenv_(integer *, char *, char *, integer *, integer *, 
954             integer *, integer *, ftnlen, ftnlen);
955     extern real slamch_(char *);
956     integer ijobvr;
957     logical wantsb;
958     extern /* Subroutine */ int cungqr_(integer *, integer *, integer *, 
959             complex *, integer *, complex *, complex *, integer *, integer *);
960     real anrmto;
961     logical wantse;
962     real bnrmto;
963     extern /* Subroutine */ int cunmqr_(char *, char *, integer *, integer *, 
964             integer *, complex *, integer *, complex *, complex *, integer *, 
965             complex *, integer *, integer *);
966     integer minwrk, maxwrk;
967     logical wantsn;
968     real smlnum;
969     logical lquery, wantsv;
970     real eps;
971     logical ilv;
972
973
974 /*  -- LAPACK driver routine (version 3.7.0) -- */
975 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
976 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
977 /*     April 2012 */
978
979
980 /*  ===================================================================== */
981
982
983 /*     Decode the input arguments */
984
985     /* Parameter adjustments */
986     a_dim1 = *lda;
987     a_offset = 1 + a_dim1 * 1;
988     a -= a_offset;
989     b_dim1 = *ldb;
990     b_offset = 1 + b_dim1 * 1;
991     b -= b_offset;
992     --alpha;
993     --beta;
994     vl_dim1 = *ldvl;
995     vl_offset = 1 + vl_dim1 * 1;
996     vl -= vl_offset;
997     vr_dim1 = *ldvr;
998     vr_offset = 1 + vr_dim1 * 1;
999     vr -= vr_offset;
1000     --lscale;
1001     --rscale;
1002     --rconde;
1003     --rcondv;
1004     --work;
1005     --rwork;
1006     --iwork;
1007     --bwork;
1008
1009     /* Function Body */
1010     if (lsame_(jobvl, "N")) {
1011         ijobvl = 1;
1012         ilvl = FALSE_;
1013     } else if (lsame_(jobvl, "V")) {
1014         ijobvl = 2;
1015         ilvl = TRUE_;
1016     } else {
1017         ijobvl = -1;
1018         ilvl = FALSE_;
1019     }
1020
1021     if (lsame_(jobvr, "N")) {
1022         ijobvr = 1;
1023         ilvr = FALSE_;
1024     } else if (lsame_(jobvr, "V")) {
1025         ijobvr = 2;
1026         ilvr = TRUE_;
1027     } else {
1028         ijobvr = -1;
1029         ilvr = FALSE_;
1030     }
1031     ilv = ilvl || ilvr;
1032
1033     noscl = lsame_(balanc, "N") || lsame_(balanc, "P");
1034     wantsn = lsame_(sense, "N");
1035     wantse = lsame_(sense, "E");
1036     wantsv = lsame_(sense, "V");
1037     wantsb = lsame_(sense, "B");
1038
1039 /*     Test the input arguments */
1040
1041     *info = 0;
1042     lquery = *lwork == -1;
1043     if (! (noscl || lsame_(balanc, "S") || lsame_(
1044             balanc, "B"))) {
1045         *info = -1;
1046     } else if (ijobvl <= 0) {
1047         *info = -2;
1048     } else if (ijobvr <= 0) {
1049         *info = -3;
1050     } else if (! (wantsn || wantse || wantsb || wantsv)) {
1051         *info = -4;
1052     } else if (*n < 0) {
1053         *info = -5;
1054     } else if (*lda < f2cmax(1,*n)) {
1055         *info = -7;
1056     } else if (*ldb < f2cmax(1,*n)) {
1057         *info = -9;
1058     } else if (*ldvl < 1 || ilvl && *ldvl < *n) {
1059         *info = -13;
1060     } else if (*ldvr < 1 || ilvr && *ldvr < *n) {
1061         *info = -15;
1062     }
1063
1064 /*     Compute workspace */
1065 /*      (Note: Comments in the code beginning "Workspace:" describe the */
1066 /*       minimal amount of workspace needed at that point in the code, */
1067 /*       as well as the preferred amount for good performance. */
1068 /*       NB refers to the optimal block size for the immediately */
1069 /*       following subroutine, as returned by ILAENV. The workspace is */
1070 /*       computed assuming ILO = 1 and IHI = N, the worst case.) */
1071
1072     if (*info == 0) {
1073         if (*n == 0) {
1074             minwrk = 1;
1075             maxwrk = 1;
1076         } else {
1077             minwrk = *n << 1;
1078             if (wantse) {
1079                 minwrk = *n << 2;
1080             } else if (wantsv || wantsb) {
1081                 minwrk = (*n << 1) * (*n + 1);
1082             }
1083             maxwrk = minwrk;
1084 /* Computing MAX */
1085             i__1 = maxwrk, i__2 = *n + *n * ilaenv_(&c__1, "CGEQRF", " ", n, &
1086                     c__1, n, &c__0, (ftnlen)6, (ftnlen)1);
1087             maxwrk = f2cmax(i__1,i__2);
1088 /* Computing MAX */
1089             i__1 = maxwrk, i__2 = *n + *n * ilaenv_(&c__1, "CUNMQR", " ", n, &
1090                     c__1, n, &c__0, (ftnlen)6, (ftnlen)1);
1091             maxwrk = f2cmax(i__1,i__2);
1092             if (ilvl) {
1093 /* Computing MAX */
1094                 i__1 = maxwrk, i__2 = *n + *n * ilaenv_(&c__1, "CUNGQR", 
1095                         " ", n, &c__1, n, &c__0, (ftnlen)6, (ftnlen)1);
1096                 maxwrk = f2cmax(i__1,i__2);
1097             }
1098         }
1099         work[1].r = (real) maxwrk, work[1].i = 0.f;
1100
1101         if (*lwork < minwrk && ! lquery) {
1102             *info = -25;
1103         }
1104     }
1105
1106     if (*info != 0) {
1107         i__1 = -(*info);
1108         xerbla_("CGGEVX", &i__1, (ftnlen)6);
1109         return 0;
1110     } else if (lquery) {
1111         return 0;
1112     }
1113
1114 /*     Quick return if possible */
1115
1116     if (*n == 0) {
1117         return 0;
1118     }
1119
1120 /*     Get machine constants */
1121
1122     eps = slamch_("P");
1123     smlnum = slamch_("S");
1124     bignum = 1.f / smlnum;
1125     slabad_(&smlnum, &bignum);
1126     smlnum = sqrt(smlnum) / eps;
1127     bignum = 1.f / smlnum;
1128
1129 /*     Scale A if f2cmax element outside range [SMLNUM,BIGNUM] */
1130
1131     anrm = clange_("M", n, n, &a[a_offset], lda, &rwork[1]);
1132     ilascl = FALSE_;
1133     if (anrm > 0.f && anrm < smlnum) {
1134         anrmto = smlnum;
1135         ilascl = TRUE_;
1136     } else if (anrm > bignum) {
1137         anrmto = bignum;
1138         ilascl = TRUE_;
1139     }
1140     if (ilascl) {
1141         clascl_("G", &c__0, &c__0, &anrm, &anrmto, n, n, &a[a_offset], lda, &
1142                 ierr);
1143     }
1144
1145 /*     Scale B if f2cmax element outside range [SMLNUM,BIGNUM] */
1146
1147     bnrm = clange_("M", n, n, &b[b_offset], ldb, &rwork[1]);
1148     ilbscl = FALSE_;
1149     if (bnrm > 0.f && bnrm < smlnum) {
1150         bnrmto = smlnum;
1151         ilbscl = TRUE_;
1152     } else if (bnrm > bignum) {
1153         bnrmto = bignum;
1154         ilbscl = TRUE_;
1155     }
1156     if (ilbscl) {
1157         clascl_("G", &c__0, &c__0, &bnrm, &bnrmto, n, n, &b[b_offset], ldb, &
1158                 ierr);
1159     }
1160
1161 /*     Permute and/or balance the matrix pair (A,B) */
1162 /*     (Real Workspace: need 6*N if BALANC = 'S' or 'B', 1 otherwise) */
1163
1164     cggbal_(balanc, n, &a[a_offset], lda, &b[b_offset], ldb, ilo, ihi, &
1165             lscale[1], &rscale[1], &rwork[1], &ierr);
1166
1167 /*     Compute ABNRM and BBNRM */
1168
1169     *abnrm = clange_("1", n, n, &a[a_offset], lda, &rwork[1]);
1170     if (ilascl) {
1171         rwork[1] = *abnrm;
1172         slascl_("G", &c__0, &c__0, &anrmto, &anrm, &c__1, &c__1, &rwork[1], &
1173                 c__1, &ierr);
1174         *abnrm = rwork[1];
1175     }
1176
1177     *bbnrm = clange_("1", n, n, &b[b_offset], ldb, &rwork[1]);
1178     if (ilbscl) {
1179         rwork[1] = *bbnrm;
1180         slascl_("G", &c__0, &c__0, &bnrmto, &bnrm, &c__1, &c__1, &rwork[1], &
1181                 c__1, &ierr);
1182         *bbnrm = rwork[1];
1183     }
1184
1185 /*     Reduce B to triangular form (QR decomposition of B) */
1186 /*     (Complex Workspace: need N, prefer N*NB ) */
1187
1188     irows = *ihi + 1 - *ilo;
1189     if (ilv || ! wantsn) {
1190         icols = *n + 1 - *ilo;
1191     } else {
1192         icols = irows;
1193     }
1194     itau = 1;
1195     iwrk = itau + irows;
1196     i__1 = *lwork + 1 - iwrk;
1197     cgeqrf_(&irows, &icols, &b[*ilo + *ilo * b_dim1], ldb, &work[itau], &work[
1198             iwrk], &i__1, &ierr);
1199
1200 /*     Apply the unitary transformation to A */
1201 /*     (Complex Workspace: need N, prefer N*NB) */
1202
1203     i__1 = *lwork + 1 - iwrk;
1204     cunmqr_("L", "C", &irows, &icols, &irows, &b[*ilo + *ilo * b_dim1], ldb, &
1205             work[itau], &a[*ilo + *ilo * a_dim1], lda, &work[iwrk], &i__1, &
1206             ierr);
1207
1208 /*     Initialize VL and/or VR */
1209 /*     (Workspace: need N, prefer N*NB) */
1210
1211     if (ilvl) {
1212         claset_("Full", n, n, &c_b1, &c_b2, &vl[vl_offset], ldvl);
1213         if (irows > 1) {
1214             i__1 = irows - 1;
1215             i__2 = irows - 1;
1216             clacpy_("L", &i__1, &i__2, &b[*ilo + 1 + *ilo * b_dim1], ldb, &vl[
1217                     *ilo + 1 + *ilo * vl_dim1], ldvl);
1218         }
1219         i__1 = *lwork + 1 - iwrk;
1220         cungqr_(&irows, &irows, &irows, &vl[*ilo + *ilo * vl_dim1], ldvl, &
1221                 work[itau], &work[iwrk], &i__1, &ierr);
1222     }
1223
1224     if (ilvr) {
1225         claset_("Full", n, n, &c_b1, &c_b2, &vr[vr_offset], ldvr);
1226     }
1227
1228 /*     Reduce to generalized Hessenberg form */
1229 /*     (Workspace: none needed) */
1230
1231     if (ilv || ! wantsn) {
1232
1233 /*        Eigenvectors requested -- work on whole matrix. */
1234
1235         cgghrd_(jobvl, jobvr, n, ilo, ihi, &a[a_offset], lda, &b[b_offset], 
1236                 ldb, &vl[vl_offset], ldvl, &vr[vr_offset], ldvr, &ierr);
1237     } else {
1238         cgghrd_("N", "N", &irows, &c__1, &irows, &a[*ilo + *ilo * a_dim1], 
1239                 lda, &b[*ilo + *ilo * b_dim1], ldb, &vl[vl_offset], ldvl, &vr[
1240                 vr_offset], ldvr, &ierr);
1241     }
1242
1243 /*     Perform QZ algorithm (Compute eigenvalues, and optionally, the */
1244 /*     Schur forms and Schur vectors) */
1245 /*     (Complex Workspace: need N) */
1246 /*     (Real Workspace: need N) */
1247
1248     iwrk = itau;
1249     if (ilv || ! wantsn) {
1250         *(unsigned char *)chtemp = 'S';
1251     } else {
1252         *(unsigned char *)chtemp = 'E';
1253     }
1254
1255     i__1 = *lwork + 1 - iwrk;
1256     chgeqz_(chtemp, jobvl, jobvr, n, ilo, ihi, &a[a_offset], lda, &b[b_offset]
1257             , ldb, &alpha[1], &beta[1], &vl[vl_offset], ldvl, &vr[vr_offset], 
1258             ldvr, &work[iwrk], &i__1, &rwork[1], &ierr);
1259     if (ierr != 0) {
1260         if (ierr > 0 && ierr <= *n) {
1261             *info = ierr;
1262         } else if (ierr > *n && ierr <= *n << 1) {
1263             *info = ierr - *n;
1264         } else {
1265             *info = *n + 1;
1266         }
1267         goto L90;
1268     }
1269
1270 /*     Compute Eigenvectors and estimate condition numbers if desired */
1271 /*     CTGEVC: (Complex Workspace: need 2*N ) */
1272 /*             (Real Workspace:    need 2*N ) */
1273 /*     CTGSNA: (Complex Workspace: need 2*N*N if SENSE='V' or 'B') */
1274 /*             (Integer Workspace: need N+2 ) */
1275
1276     if (ilv || ! wantsn) {
1277         if (ilv) {
1278             if (ilvl) {
1279                 if (ilvr) {
1280                     *(unsigned char *)chtemp = 'B';
1281                 } else {
1282                     *(unsigned char *)chtemp = 'L';
1283                 }
1284             } else {
1285                 *(unsigned char *)chtemp = 'R';
1286             }
1287
1288             ctgevc_(chtemp, "B", ldumma, n, &a[a_offset], lda, &b[b_offset], 
1289                     ldb, &vl[vl_offset], ldvl, &vr[vr_offset], ldvr, n, &in, &
1290                     work[iwrk], &rwork[1], &ierr);
1291             if (ierr != 0) {
1292                 *info = *n + 2;
1293                 goto L90;
1294             }
1295         }
1296
1297         if (! wantsn) {
1298
1299 /*           compute eigenvectors (STGEVC) and estimate condition */
1300 /*           numbers (STGSNA). Note that the definition of the condition */
1301 /*           number is not invariant under transformation (u,v) to */
1302 /*           (Q*u, Z*v), where (u,v) are eigenvectors of the generalized */
1303 /*           Schur form (S,T), Q and Z are orthogonal matrices. In order */
1304 /*           to avoid using extra 2*N*N workspace, we have to */
1305 /*           re-calculate eigenvectors and estimate the condition numbers */
1306 /*           one at a time. */
1307
1308             i__1 = *n;
1309             for (i__ = 1; i__ <= i__1; ++i__) {
1310
1311                 i__2 = *n;
1312                 for (j = 1; j <= i__2; ++j) {
1313                     bwork[j] = FALSE_;
1314 /* L10: */
1315                 }
1316                 bwork[i__] = TRUE_;
1317
1318                 iwrk = *n + 1;
1319                 iwrk1 = iwrk + *n;
1320
1321                 if (wantse || wantsb) {
1322                     ctgevc_("B", "S", &bwork[1], n, &a[a_offset], lda, &b[
1323                             b_offset], ldb, &work[1], n, &work[iwrk], n, &
1324                             c__1, &m, &work[iwrk1], &rwork[1], &ierr);
1325                     if (ierr != 0) {
1326                         *info = *n + 2;
1327                         goto L90;
1328                     }
1329                 }
1330
1331                 i__2 = *lwork - iwrk1 + 1;
1332                 ctgsna_(sense, "S", &bwork[1], n, &a[a_offset], lda, &b[
1333                         b_offset], ldb, &work[1], n, &work[iwrk], n, &rconde[
1334                         i__], &rcondv[i__], &c__1, &m, &work[iwrk1], &i__2, &
1335                         iwork[1], &ierr);
1336
1337 /* L20: */
1338             }
1339         }
1340     }
1341
1342 /*     Undo balancing on VL and VR and normalization */
1343 /*     (Workspace: none needed) */
1344
1345     if (ilvl) {
1346         cggbak_(balanc, "L", n, ilo, ihi, &lscale[1], &rscale[1], n, &vl[
1347                 vl_offset], ldvl, &ierr);
1348
1349         i__1 = *n;
1350         for (jc = 1; jc <= i__1; ++jc) {
1351             temp = 0.f;
1352             i__2 = *n;
1353             for (jr = 1; jr <= i__2; ++jr) {
1354 /* Computing MAX */
1355                 i__3 = jr + jc * vl_dim1;
1356                 r__3 = temp, r__4 = (r__1 = vl[i__3].r, abs(r__1)) + (r__2 = 
1357                         r_imag(&vl[jr + jc * vl_dim1]), abs(r__2));
1358                 temp = f2cmax(r__3,r__4);
1359 /* L30: */
1360             }
1361             if (temp < smlnum) {
1362                 goto L50;
1363             }
1364             temp = 1.f / temp;
1365             i__2 = *n;
1366             for (jr = 1; jr <= i__2; ++jr) {
1367                 i__3 = jr + jc * vl_dim1;
1368                 i__4 = jr + jc * vl_dim1;
1369                 q__1.r = temp * vl[i__4].r, q__1.i = temp * vl[i__4].i;
1370                 vl[i__3].r = q__1.r, vl[i__3].i = q__1.i;
1371 /* L40: */
1372             }
1373 L50:
1374             ;
1375         }
1376     }
1377
1378     if (ilvr) {
1379         cggbak_(balanc, "R", n, ilo, ihi, &lscale[1], &rscale[1], n, &vr[
1380                 vr_offset], ldvr, &ierr);
1381         i__1 = *n;
1382         for (jc = 1; jc <= i__1; ++jc) {
1383             temp = 0.f;
1384             i__2 = *n;
1385             for (jr = 1; jr <= i__2; ++jr) {
1386 /* Computing MAX */
1387                 i__3 = jr + jc * vr_dim1;
1388                 r__3 = temp, r__4 = (r__1 = vr[i__3].r, abs(r__1)) + (r__2 = 
1389                         r_imag(&vr[jr + jc * vr_dim1]), abs(r__2));
1390                 temp = f2cmax(r__3,r__4);
1391 /* L60: */
1392             }
1393             if (temp < smlnum) {
1394                 goto L80;
1395             }
1396             temp = 1.f / temp;
1397             i__2 = *n;
1398             for (jr = 1; jr <= i__2; ++jr) {
1399                 i__3 = jr + jc * vr_dim1;
1400                 i__4 = jr + jc * vr_dim1;
1401                 q__1.r = temp * vr[i__4].r, q__1.i = temp * vr[i__4].i;
1402                 vr[i__3].r = q__1.r, vr[i__3].i = q__1.i;
1403 /* L70: */
1404             }
1405 L80:
1406             ;
1407         }
1408     }
1409
1410 /*     Undo scaling if necessary */
1411
1412 L90:
1413
1414     if (ilascl) {
1415         clascl_("G", &c__0, &c__0, &anrmto, &anrm, n, &c__1, &alpha[1], n, &
1416                 ierr);
1417     }
1418
1419     if (ilbscl) {
1420         clascl_("G", &c__0, &c__0, &bnrmto, &bnrm, n, &c__1, &beta[1], n, &
1421                 ierr);
1422     }
1423
1424     work[1].r = (real) maxwrk, work[1].i = 0.f;
1425     return 0;
1426
1427 /*     End of CGGEVX */
1428
1429 } /* cggevx_ */
1430