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