C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / zgges3.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 doublecomplex c_b1 = {0.,0.};
516 static doublecomplex c_b2 = {1.,0.};
517 static integer c_n1 = -1;
518 static integer c__1 = 1;
519 static integer c__0 = 0;
520
521 /* > \brief <b> ZGGES3 computes the eigenvalues, the Schur form, and, optionally, the matrix of Schur vectors 
522 for GE matrices (blocked algorithm)</b> */
523
524 /*  =========== DOCUMENTATION =========== */
525
526 /* Online html documentation available at */
527 /*            http://www.netlib.org/lapack/explore-html/ */
528
529 /* > \htmlonly */
530 /* > Download ZGGES3 + dependencies */
531 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgges3.
532 f"> */
533 /* > [TGZ]</a> */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgges3.
535 f"> */
536 /* > [ZIP]</a> */
537 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgges3.
538 f"> */
539 /* > [TXT]</a> */
540 /* > \endhtmlonly */
541
542 /*  Definition: */
543 /*  =========== */
544
545 /*       SUBROUTINE ZGGES3( JOBVSL, JOBVSR, SORT, SELCTG, N, A, LDA, B, */
546 /*      $                   LDB, SDIM, ALPHA, BETA, VSL, LDVSL, VSR, LDVSR, */
547 /*      $                   WORK, LWORK, RWORK, BWORK, INFO ) */
548
549 /*       CHARACTER          JOBVSL, JOBVSR, SORT */
550 /*       INTEGER            INFO, LDA, LDB, LDVSL, LDVSR, LWORK, N, SDIM */
551 /*       LOGICAL            BWORK( * ) */
552 /*       DOUBLE PRECISION   RWORK( * ) */
553 /*       COMPLEX*16         A( LDA, * ), ALPHA( * ), B( LDB, * ), */
554 /*      $                   BETA( * ), VSL( LDVSL, * ), VSR( LDVSR, * ), */
555 /*      $                   WORK( * ) */
556 /*       LOGICAL            SELCTG */
557 /*       EXTERNAL           SELCTG */
558
559
560 /* > \par Purpose: */
561 /*  ============= */
562 /* > */
563 /* > \verbatim */
564 /* > */
565 /* > ZGGES3 computes for a pair of N-by-N complex nonsymmetric matrices */
566 /* > (A,B), the generalized eigenvalues, the generalized complex Schur */
567 /* > form (S, T), and optionally left and/or right Schur vectors (VSL */
568 /* > and VSR). This gives the generalized Schur factorization */
569 /* > */
570 /* >         (A,B) = ( (VSL)*S*(VSR)**H, (VSL)*T*(VSR)**H ) */
571 /* > */
572 /* > where (VSR)**H is the conjugate-transpose of VSR. */
573 /* > */
574 /* > Optionally, it also orders the eigenvalues so that a selected cluster */
575 /* > of eigenvalues appears in the leading diagonal blocks of the upper */
576 /* > triangular matrix S and the upper triangular matrix T. The leading */
577 /* > columns of VSL and VSR then form an unitary basis for the */
578 /* > corresponding left and right eigenspaces (deflating subspaces). */
579 /* > */
580 /* > (If only the generalized eigenvalues are needed, use the driver */
581 /* > ZGGEV instead, which is faster.) */
582 /* > */
583 /* > A generalized eigenvalue for a pair of matrices (A,B) is a scalar w */
584 /* > or a ratio alpha/beta = w, such that  A - w*B is singular.  It is */
585 /* > usually represented as the pair (alpha,beta), as there is a */
586 /* > reasonable interpretation for beta=0, and even for both being zero. */
587 /* > */
588 /* > A pair of matrices (S,T) is in generalized complex Schur form if S */
589 /* > and T are upper triangular and, in addition, the diagonal elements */
590 /* > of T are non-negative real numbers. */
591 /* > \endverbatim */
592
593 /*  Arguments: */
594 /*  ========== */
595
596 /* > \param[in] JOBVSL */
597 /* > \verbatim */
598 /* >          JOBVSL is CHARACTER*1 */
599 /* >          = 'N':  do not compute the left Schur vectors; */
600 /* >          = 'V':  compute the left Schur vectors. */
601 /* > \endverbatim */
602 /* > */
603 /* > \param[in] JOBVSR */
604 /* > \verbatim */
605 /* >          JOBVSR is CHARACTER*1 */
606 /* >          = 'N':  do not compute the right Schur vectors; */
607 /* >          = 'V':  compute the right Schur vectors. */
608 /* > \endverbatim */
609 /* > */
610 /* > \param[in] SORT */
611 /* > \verbatim */
612 /* >          SORT is CHARACTER*1 */
613 /* >          Specifies whether or not to order the eigenvalues on the */
614 /* >          diagonal of the generalized Schur form. */
615 /* >          = 'N':  Eigenvalues are not ordered; */
616 /* >          = 'S':  Eigenvalues are ordered (see SELCTG). */
617 /* > \endverbatim */
618 /* > */
619 /* > \param[in] SELCTG */
620 /* > \verbatim */
621 /* >          SELCTG is a LOGICAL FUNCTION of two COMPLEX*16 arguments */
622 /* >          SELCTG must be declared EXTERNAL in the calling subroutine. */
623 /* >          If SORT = 'N', SELCTG is not referenced. */
624 /* >          If SORT = 'S', SELCTG is used to select eigenvalues to sort */
625 /* >          to the top left of the Schur form. */
626 /* >          An eigenvalue ALPHA(j)/BETA(j) is selected if */
627 /* >          SELCTG(ALPHA(j),BETA(j)) is true. */
628 /* > */
629 /* >          Note that a selected complex eigenvalue may no longer satisfy */
630 /* >          SELCTG(ALPHA(j),BETA(j)) = .TRUE. after ordering, since */
631 /* >          ordering may change the value of complex eigenvalues */
632 /* >          (especially if the eigenvalue is ill-conditioned), in this */
633 /* >          case INFO is set to N+2 (See INFO below). */
634 /* > \endverbatim */
635 /* > */
636 /* > \param[in] N */
637 /* > \verbatim */
638 /* >          N is INTEGER */
639 /* >          The order of the matrices A, B, VSL, and VSR.  N >= 0. */
640 /* > \endverbatim */
641 /* > */
642 /* > \param[in,out] A */
643 /* > \verbatim */
644 /* >          A is COMPLEX*16 array, dimension (LDA, N) */
645 /* >          On entry, the first of the pair of matrices. */
646 /* >          On exit, A has been overwritten by its generalized Schur */
647 /* >          form S. */
648 /* > \endverbatim */
649 /* > */
650 /* > \param[in] LDA */
651 /* > \verbatim */
652 /* >          LDA is INTEGER */
653 /* >          The leading dimension of A.  LDA >= f2cmax(1,N). */
654 /* > \endverbatim */
655 /* > */
656 /* > \param[in,out] B */
657 /* > \verbatim */
658 /* >          B is COMPLEX*16 array, dimension (LDB, N) */
659 /* >          On entry, the second of the pair of matrices. */
660 /* >          On exit, B has been overwritten by its generalized Schur */
661 /* >          form T. */
662 /* > \endverbatim */
663 /* > */
664 /* > \param[in] LDB */
665 /* > \verbatim */
666 /* >          LDB is INTEGER */
667 /* >          The leading dimension of B.  LDB >= f2cmax(1,N). */
668 /* > \endverbatim */
669 /* > */
670 /* > \param[out] SDIM */
671 /* > \verbatim */
672 /* >          SDIM is INTEGER */
673 /* >          If SORT = 'N', SDIM = 0. */
674 /* >          If SORT = 'S', SDIM = number of eigenvalues (after sorting) */
675 /* >          for which SELCTG is true. */
676 /* > \endverbatim */
677 /* > */
678 /* > \param[out] ALPHA */
679 /* > \verbatim */
680 /* >          ALPHA is COMPLEX*16 array, dimension (N) */
681 /* > \endverbatim */
682 /* > */
683 /* > \param[out] BETA */
684 /* > \verbatim */
685 /* >          BETA is COMPLEX*16 array, dimension (N) */
686 /* >          On exit,  ALPHA(j)/BETA(j), j=1,...,N, will be the */
687 /* >          generalized eigenvalues.  ALPHA(j), j=1,...,N  and  BETA(j), */
688 /* >          j=1,...,N  are the diagonals of the complex Schur form (A,B) */
689 /* >          output by ZGGES3. The  BETA(j) will be non-negative real. */
690 /* > */
691 /* >          Note: the quotients ALPHA(j)/BETA(j) may easily over- or */
692 /* >          underflow, and BETA(j) may even be zero.  Thus, the user */
693 /* >          should avoid naively computing the ratio alpha/beta. */
694 /* >          However, ALPHA will be always less than and usually */
695 /* >          comparable with norm(A) in magnitude, and BETA always less */
696 /* >          than and usually comparable with norm(B). */
697 /* > \endverbatim */
698 /* > */
699 /* > \param[out] VSL */
700 /* > \verbatim */
701 /* >          VSL is COMPLEX*16 array, dimension (LDVSL,N) */
702 /* >          If JOBVSL = 'V', VSL will contain the left Schur vectors. */
703 /* >          Not referenced if JOBVSL = 'N'. */
704 /* > \endverbatim */
705 /* > */
706 /* > \param[in] LDVSL */
707 /* > \verbatim */
708 /* >          LDVSL is INTEGER */
709 /* >          The leading dimension of the matrix VSL. LDVSL >= 1, and */
710 /* >          if JOBVSL = 'V', LDVSL >= N. */
711 /* > \endverbatim */
712 /* > */
713 /* > \param[out] VSR */
714 /* > \verbatim */
715 /* >          VSR is COMPLEX*16 array, dimension (LDVSR,N) */
716 /* >          If JOBVSR = 'V', VSR will contain the right Schur vectors. */
717 /* >          Not referenced if JOBVSR = 'N'. */
718 /* > \endverbatim */
719 /* > */
720 /* > \param[in] LDVSR */
721 /* > \verbatim */
722 /* >          LDVSR is INTEGER */
723 /* >          The leading dimension of the matrix VSR. LDVSR >= 1, and */
724 /* >          if JOBVSR = 'V', LDVSR >= N. */
725 /* > \endverbatim */
726 /* > */
727 /* > \param[out] WORK */
728 /* > \verbatim */
729 /* >          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) */
730 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
731 /* > \endverbatim */
732 /* > */
733 /* > \param[in] LWORK */
734 /* > \verbatim */
735 /* >          LWORK is INTEGER */
736 /* >          The dimension of the array WORK. */
737 /* > */
738 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
739 /* >          only calculates the optimal size of the WORK array, returns */
740 /* >          this value as the first entry of the WORK array, and no error */
741 /* >          message related to LWORK is issued by XERBLA. */
742 /* > \endverbatim */
743 /* > */
744 /* > \param[out] RWORK */
745 /* > \verbatim */
746 /* >          RWORK is DOUBLE PRECISION array, dimension (8*N) */
747 /* > \endverbatim */
748 /* > */
749 /* > \param[out] BWORK */
750 /* > \verbatim */
751 /* >          BWORK is LOGICAL array, dimension (N) */
752 /* >          Not referenced if SORT = 'N'. */
753 /* > \endverbatim */
754 /* > */
755 /* > \param[out] INFO */
756 /* > \verbatim */
757 /* >          INFO is INTEGER */
758 /* >          = 0:  successful exit */
759 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value. */
760 /* >          =1,...,N: */
761 /* >                The QZ iteration failed.  (A,B) are not in Schur */
762 /* >                form, but ALPHA(j) and BETA(j) should be correct for */
763 /* >                j=INFO+1,...,N. */
764 /* >          > N:  =N+1: other than QZ iteration failed in ZHGEQZ */
765 /* >                =N+2: after reordering, roundoff changed values of */
766 /* >                      some complex eigenvalues so that leading */
767 /* >                      eigenvalues in the Generalized Schur form no */
768 /* >                      longer satisfy SELCTG=.TRUE.  This could also */
769 /* >                      be caused due to scaling. */
770 /* >                =N+3: reordering failed in ZTGSEN. */
771 /* > \endverbatim */
772
773 /*  Authors: */
774 /*  ======== */
775
776 /* > \author Univ. of Tennessee */
777 /* > \author Univ. of California Berkeley */
778 /* > \author Univ. of Colorado Denver */
779 /* > \author NAG Ltd. */
780
781 /* > \date January 2015 */
782
783 /* > \ingroup complex16GEeigen */
784
785 /*  ===================================================================== */
786 /* Subroutine */ int zgges3_(char *jobvsl, char *jobvsr, char *sort, L_fp 
787         selctg, integer *n, doublecomplex *a, integer *lda, doublecomplex *b, 
788         integer *ldb, integer *sdim, doublecomplex *alpha, doublecomplex *
789         beta, doublecomplex *vsl, integer *ldvsl, doublecomplex *vsr, integer 
790         *ldvsr, doublecomplex *work, integer *lwork, doublereal *rwork, 
791         logical *bwork, integer *info)
792 {
793     /* System generated locals */
794     integer a_dim1, a_offset, b_dim1, b_offset, vsl_dim1, vsl_offset, 
795             vsr_dim1, vsr_offset, i__1, i__2;
796     doublecomplex z__1;
797
798     /* Local variables */
799     doublereal anrm, bnrm;
800     integer idum[1], ierr, itau, iwrk;
801     doublereal pvsl, pvsr;
802     integer i__;
803     extern logical lsame_(char *, char *);
804     integer ileft, icols;
805     logical cursl, ilvsl, ilvsr;
806     integer irwrk, irows;
807     extern /* Subroutine */ int zgghd3_(char *, char *, integer *, integer *, 
808             integer *, doublecomplex *, integer *, doublecomplex *, integer *,
809              doublecomplex *, integer *, doublecomplex *, integer *, 
810             doublecomplex *, integer *, integer *), dlabad_(
811             doublereal *, doublereal *);
812     extern doublereal dlamch_(char *);
813     extern /* Subroutine */ int zggbak_(char *, char *, integer *, integer *, 
814             integer *, doublereal *, doublereal *, integer *, doublecomplex *,
815              integer *, integer *), zggbal_(char *, integer *,
816              doublecomplex *, integer *, doublecomplex *, integer *, integer *
817             , integer *, doublereal *, doublereal *, doublereal *, integer *);
818     logical ilascl, ilbscl;
819     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
820     extern doublereal zlange_(char *, integer *, integer *, doublecomplex *, 
821             integer *, doublereal *);
822     doublereal bignum;
823     integer ijobvl, iright;
824     extern /* Subroutine */ int zlascl_(char *, integer *, integer *, 
825             doublereal *, doublereal *, integer *, integer *, doublecomplex *,
826              integer *, integer *);
827     integer ijobvr;
828     extern /* Subroutine */ int zgeqrf_(integer *, integer *, doublecomplex *,
829              integer *, doublecomplex *, doublecomplex *, integer *, integer *
830             );
831     doublereal anrmto, bnrmto;
832     logical lastsl;
833     extern /* Subroutine */ int zlacpy_(char *, integer *, integer *, 
834             doublecomplex *, integer *, doublecomplex *, integer *), 
835             zlaset_(char *, integer *, integer *, doublecomplex *, 
836             doublecomplex *, doublecomplex *, integer *), zhgeqz_(
837             char *, char *, char *, integer *, integer *, integer *, 
838             doublecomplex *, integer *, doublecomplex *, integer *, 
839             doublecomplex *, doublecomplex *, doublecomplex *, integer *, 
840             doublecomplex *, integer *, doublecomplex *, integer *, 
841             doublereal *, integer *), ztgsen_(integer 
842             *, logical *, logical *, logical *, integer *, doublecomplex *, 
843             integer *, doublecomplex *, integer *, doublecomplex *, 
844             doublecomplex *, doublecomplex *, integer *, doublecomplex *, 
845             integer *, integer *, doublereal *, doublereal *, doublereal *, 
846             doublecomplex *, integer *, integer *, integer *, integer *);
847     doublereal smlnum;
848     logical wantst, lquery;
849     integer lwkopt;
850     extern /* Subroutine */ int zungqr_(integer *, integer *, integer *, 
851             doublecomplex *, integer *, doublecomplex *, doublecomplex *, 
852             integer *, integer *), zunmqr_(char *, char *, integer *, integer 
853             *, integer *, doublecomplex *, integer *, doublecomplex *, 
854             doublecomplex *, integer *, doublecomplex *, integer *, integer *);
855     doublereal dif[2];
856     integer ihi, ilo;
857     doublereal eps;
858
859
860 /*  -- LAPACK driver routine (version 3.6.1) -- */
861 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
862 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
863 /*     January 2015 */
864
865
866 /*  ===================================================================== */
867
868
869 /*     Decode the input arguments */
870
871     /* Parameter adjustments */
872     a_dim1 = *lda;
873     a_offset = 1 + a_dim1 * 1;
874     a -= a_offset;
875     b_dim1 = *ldb;
876     b_offset = 1 + b_dim1 * 1;
877     b -= b_offset;
878     --alpha;
879     --beta;
880     vsl_dim1 = *ldvsl;
881     vsl_offset = 1 + vsl_dim1 * 1;
882     vsl -= vsl_offset;
883     vsr_dim1 = *ldvsr;
884     vsr_offset = 1 + vsr_dim1 * 1;
885     vsr -= vsr_offset;
886     --work;
887     --rwork;
888     --bwork;
889
890     /* Function Body */
891     if (lsame_(jobvsl, "N")) {
892         ijobvl = 1;
893         ilvsl = FALSE_;
894     } else if (lsame_(jobvsl, "V")) {
895         ijobvl = 2;
896         ilvsl = TRUE_;
897     } else {
898         ijobvl = -1;
899         ilvsl = FALSE_;
900     }
901
902     if (lsame_(jobvsr, "N")) {
903         ijobvr = 1;
904         ilvsr = FALSE_;
905     } else if (lsame_(jobvsr, "V")) {
906         ijobvr = 2;
907         ilvsr = TRUE_;
908     } else {
909         ijobvr = -1;
910         ilvsr = FALSE_;
911     }
912
913     wantst = lsame_(sort, "S");
914
915 /*     Test the input arguments */
916
917     *info = 0;
918     lquery = *lwork == -1;
919     if (ijobvl <= 0) {
920         *info = -1;
921     } else if (ijobvr <= 0) {
922         *info = -2;
923     } else if (! wantst && ! lsame_(sort, "N")) {
924         *info = -3;
925     } else if (*n < 0) {
926         *info = -5;
927     } else if (*lda < f2cmax(1,*n)) {
928         *info = -7;
929     } else if (*ldb < f2cmax(1,*n)) {
930         *info = -9;
931     } else if (*ldvsl < 1 || ilvsl && *ldvsl < *n) {
932         *info = -14;
933     } else if (*ldvsr < 1 || ilvsr && *ldvsr < *n) {
934         *info = -16;
935     } else /* if(complicated condition) */ {
936 /* Computing MAX */
937         i__1 = 1, i__2 = *n << 1;
938         if (*lwork < f2cmax(i__1,i__2) && ! lquery) {
939             *info = -18;
940         }
941     }
942
943 /*     Compute workspace */
944
945     if (*info == 0) {
946         zgeqrf_(n, n, &b[b_offset], ldb, &work[1], &work[1], &c_n1, &ierr);
947 /* Computing MAX */
948         i__1 = 1, i__2 = *n + (integer) work[1].r;
949         lwkopt = f2cmax(i__1,i__2);
950         zunmqr_("L", "C", n, n, n, &b[b_offset], ldb, &work[1], &a[a_offset], 
951                 lda, &work[1], &c_n1, &ierr);
952 /* Computing MAX */
953         i__1 = lwkopt, i__2 = *n + (integer) work[1].r;
954         lwkopt = f2cmax(i__1,i__2);
955         if (ilvsl) {
956             zungqr_(n, n, n, &vsl[vsl_offset], ldvsl, &work[1], &work[1], &
957                     c_n1, &ierr);
958 /* Computing MAX */
959             i__1 = lwkopt, i__2 = *n + (integer) work[1].r;
960             lwkopt = f2cmax(i__1,i__2);
961         }
962         zgghd3_(jobvsl, jobvsr, n, &c__1, n, &a[a_offset], lda, &b[b_offset], 
963                 ldb, &vsl[vsl_offset], ldvsl, &vsr[vsr_offset], ldvsr, &work[
964                 1], &c_n1, &ierr);
965 /* Computing MAX */
966         i__1 = lwkopt, i__2 = *n + (integer) work[1].r;
967         lwkopt = f2cmax(i__1,i__2);
968         zhgeqz_("S", jobvsl, jobvsr, n, &c__1, n, &a[a_offset], lda, &b[
969                 b_offset], ldb, &alpha[1], &beta[1], &vsl[vsl_offset], ldvsl, 
970                 &vsr[vsr_offset], ldvsr, &work[1], &c_n1, &rwork[1], &ierr);
971 /* Computing MAX */
972         i__1 = lwkopt, i__2 = (integer) work[1].r;
973         lwkopt = f2cmax(i__1,i__2);
974         if (wantst) {
975             ztgsen_(&c__0, &ilvsl, &ilvsr, &bwork[1], n, &a[a_offset], lda, &
976                     b[b_offset], ldb, &alpha[1], &beta[1], &vsl[vsl_offset], 
977                     ldvsl, &vsr[vsr_offset], ldvsr, sdim, &pvsl, &pvsr, dif, &
978                     work[1], &c_n1, idum, &c__1, &ierr);
979 /* Computing MAX */
980             i__1 = lwkopt, i__2 = (integer) work[1].r;
981             lwkopt = f2cmax(i__1,i__2);
982         }
983         z__1.r = (doublereal) lwkopt, z__1.i = 0.;
984         work[1].r = z__1.r, work[1].i = z__1.i;
985     }
986
987     if (*info != 0) {
988         i__1 = -(*info);
989         xerbla_("ZGGES3 ", &i__1, (ftnlen)6);
990         return 0;
991     } else if (lquery) {
992         return 0;
993     }
994
995 /*     Quick return if possible */
996
997     if (*n == 0) {
998         *sdim = 0;
999         return 0;
1000     }
1001
1002 /*     Get machine constants */
1003
1004     eps = dlamch_("P");
1005     smlnum = dlamch_("S");
1006     bignum = 1. / smlnum;
1007     dlabad_(&smlnum, &bignum);
1008     smlnum = sqrt(smlnum) / eps;
1009     bignum = 1. / smlnum;
1010
1011 /*     Scale A if f2cmax element outside range [SMLNUM,BIGNUM] */
1012
1013     anrm = zlange_("M", n, n, &a[a_offset], lda, &rwork[1]);
1014     ilascl = FALSE_;
1015     if (anrm > 0. && anrm < smlnum) {
1016         anrmto = smlnum;
1017         ilascl = TRUE_;
1018     } else if (anrm > bignum) {
1019         anrmto = bignum;
1020         ilascl = TRUE_;
1021     }
1022
1023     if (ilascl) {
1024         zlascl_("G", &c__0, &c__0, &anrm, &anrmto, n, n, &a[a_offset], lda, &
1025                 ierr);
1026     }
1027
1028 /*     Scale B if f2cmax element outside range [SMLNUM,BIGNUM] */
1029
1030     bnrm = zlange_("M", n, n, &b[b_offset], ldb, &rwork[1]);
1031     ilbscl = FALSE_;
1032     if (bnrm > 0. && bnrm < smlnum) {
1033         bnrmto = smlnum;
1034         ilbscl = TRUE_;
1035     } else if (bnrm > bignum) {
1036         bnrmto = bignum;
1037         ilbscl = TRUE_;
1038     }
1039
1040     if (ilbscl) {
1041         zlascl_("G", &c__0, &c__0, &bnrm, &bnrmto, n, n, &b[b_offset], ldb, &
1042                 ierr);
1043     }
1044
1045 /*     Permute the matrix to make it more nearly triangular */
1046
1047     ileft = 1;
1048     iright = *n + 1;
1049     irwrk = iright + *n;
1050     zggbal_("P", n, &a[a_offset], lda, &b[b_offset], ldb, &ilo, &ihi, &rwork[
1051             ileft], &rwork[iright], &rwork[irwrk], &ierr);
1052
1053 /*     Reduce B to triangular form (QR decomposition of B) */
1054
1055     irows = ihi + 1 - ilo;
1056     icols = *n + 1 - ilo;
1057     itau = 1;
1058     iwrk = itau + irows;
1059     i__1 = *lwork + 1 - iwrk;
1060     zgeqrf_(&irows, &icols, &b[ilo + ilo * b_dim1], ldb, &work[itau], &work[
1061             iwrk], &i__1, &ierr);
1062
1063 /*     Apply the orthogonal transformation to matrix A */
1064
1065     i__1 = *lwork + 1 - iwrk;
1066     zunmqr_("L", "C", &irows, &icols, &irows, &b[ilo + ilo * b_dim1], ldb, &
1067             work[itau], &a[ilo + ilo * a_dim1], lda, &work[iwrk], &i__1, &
1068             ierr);
1069
1070 /*     Initialize VSL */
1071
1072     if (ilvsl) {
1073         zlaset_("Full", n, n, &c_b1, &c_b2, &vsl[vsl_offset], ldvsl);
1074         if (irows > 1) {
1075             i__1 = irows - 1;
1076             i__2 = irows - 1;
1077             zlacpy_("L", &i__1, &i__2, &b[ilo + 1 + ilo * b_dim1], ldb, &vsl[
1078                     ilo + 1 + ilo * vsl_dim1], ldvsl);
1079         }
1080         i__1 = *lwork + 1 - iwrk;
1081         zungqr_(&irows, &irows, &irows, &vsl[ilo + ilo * vsl_dim1], ldvsl, &
1082                 work[itau], &work[iwrk], &i__1, &ierr);
1083     }
1084
1085 /*     Initialize VSR */
1086
1087     if (ilvsr) {
1088         zlaset_("Full", n, n, &c_b1, &c_b2, &vsr[vsr_offset], ldvsr);
1089     }
1090
1091 /*     Reduce to generalized Hessenberg form */
1092
1093     i__1 = *lwork + 1 - iwrk;
1094     zgghd3_(jobvsl, jobvsr, n, &ilo, &ihi, &a[a_offset], lda, &b[b_offset], 
1095             ldb, &vsl[vsl_offset], ldvsl, &vsr[vsr_offset], ldvsr, &work[iwrk]
1096             , &i__1, &ierr);
1097
1098     *sdim = 0;
1099
1100 /*     Perform QZ algorithm, computing Schur vectors if desired */
1101
1102     iwrk = itau;
1103     i__1 = *lwork + 1 - iwrk;
1104     zhgeqz_("S", jobvsl, jobvsr, n, &ilo, &ihi, &a[a_offset], lda, &b[
1105             b_offset], ldb, &alpha[1], &beta[1], &vsl[vsl_offset], ldvsl, &
1106             vsr[vsr_offset], ldvsr, &work[iwrk], &i__1, &rwork[irwrk], &ierr);
1107     if (ierr != 0) {
1108         if (ierr > 0 && ierr <= *n) {
1109             *info = ierr;
1110         } else if (ierr > *n && ierr <= *n << 1) {
1111             *info = ierr - *n;
1112         } else {
1113             *info = *n + 1;
1114         }
1115         goto L30;
1116     }
1117
1118 /*     Sort eigenvalues ALPHA/BETA if desired */
1119
1120     if (wantst) {
1121
1122 /*        Undo scaling on eigenvalues before selecting */
1123
1124         if (ilascl) {
1125             zlascl_("G", &c__0, &c__0, &anrm, &anrmto, n, &c__1, &alpha[1], n,
1126                      &ierr);
1127         }
1128         if (ilbscl) {
1129             zlascl_("G", &c__0, &c__0, &bnrm, &bnrmto, n, &c__1, &beta[1], n, 
1130                     &ierr);
1131         }
1132
1133 /*        Select eigenvalues */
1134
1135         i__1 = *n;
1136         for (i__ = 1; i__ <= i__1; ++i__) {
1137             bwork[i__] = (*selctg)(&alpha[i__], &beta[i__]);
1138 /* L10: */
1139         }
1140
1141         i__1 = *lwork - iwrk + 1;
1142         ztgsen_(&c__0, &ilvsl, &ilvsr, &bwork[1], n, &a[a_offset], lda, &b[
1143                 b_offset], ldb, &alpha[1], &beta[1], &vsl[vsl_offset], ldvsl, 
1144                 &vsr[vsr_offset], ldvsr, sdim, &pvsl, &pvsr, dif, &work[iwrk],
1145                  &i__1, idum, &c__1, &ierr);
1146         if (ierr == 1) {
1147             *info = *n + 3;
1148         }
1149
1150     }
1151
1152 /*     Apply back-permutation to VSL and VSR */
1153
1154     if (ilvsl) {
1155         zggbak_("P", "L", n, &ilo, &ihi, &rwork[ileft], &rwork[iright], n, &
1156                 vsl[vsl_offset], ldvsl, &ierr);
1157     }
1158     if (ilvsr) {
1159         zggbak_("P", "R", n, &ilo, &ihi, &rwork[ileft], &rwork[iright], n, &
1160                 vsr[vsr_offset], ldvsr, &ierr);
1161     }
1162
1163 /*     Undo scaling */
1164
1165     if (ilascl) {
1166         zlascl_("U", &c__0, &c__0, &anrmto, &anrm, n, n, &a[a_offset], lda, &
1167                 ierr);
1168         zlascl_("G", &c__0, &c__0, &anrmto, &anrm, n, &c__1, &alpha[1], n, &
1169                 ierr);
1170     }
1171
1172     if (ilbscl) {
1173         zlascl_("U", &c__0, &c__0, &bnrmto, &bnrm, n, n, &b[b_offset], ldb, &
1174                 ierr);
1175         zlascl_("G", &c__0, &c__0, &bnrmto, &bnrm, n, &c__1, &beta[1], n, &
1176                 ierr);
1177     }
1178
1179     if (wantst) {
1180
1181 /*        Check if reordering is correct */
1182
1183         lastsl = TRUE_;
1184         *sdim = 0;
1185         i__1 = *n;
1186         for (i__ = 1; i__ <= i__1; ++i__) {
1187             cursl = (*selctg)(&alpha[i__], &beta[i__]);
1188             if (cursl) {
1189                 ++(*sdim);
1190             }
1191             if (cursl && ! lastsl) {
1192                 *info = *n + 2;
1193             }
1194             lastsl = cursl;
1195 /* L20: */
1196         }
1197
1198     }
1199
1200 L30:
1201
1202     z__1.r = (doublereal) lwkopt, z__1.i = 0.;
1203     work[1].r = z__1.r, work[1].i = z__1.i;
1204
1205     return 0;
1206
1207 /*     End of ZGGES3 */
1208
1209 } /* zgges3_ */
1210