C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / zgeev.c
1 #include <math.h>
2 #include <stdlib.h>
3 #include <string.h>
4 #include <stdio.h>
5 #include <complex.h>
6 #ifdef complex
7 #undef complex
8 #endif
9 #ifdef I
10 #undef I
11 #endif
12
13 #if defined(_WIN64)
14 typedef long long BLASLONG;
15 typedef unsigned long long BLASULONG;
16 #else
17 typedef long BLASLONG;
18 typedef unsigned long BLASULONG;
19 #endif
20
21 #ifdef LAPACK_ILP64
22 typedef BLASLONG blasint;
23 #if defined(_WIN64)
24 #define blasabs(x) llabs(x)
25 #else
26 #define blasabs(x) labs(x)
27 #endif
28 #else
29 typedef int blasint;
30 #define blasabs(x) abs(x)
31 #endif
32
33 typedef blasint integer;
34
35 typedef unsigned int uinteger;
36 typedef char *address;
37 typedef short int shortint;
38 typedef float real;
39 typedef double doublereal;
40 typedef struct { real r, i; } complex;
41 typedef struct { doublereal r, i; } doublecomplex;
42 #ifdef _MSC_VER
43 static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
44 static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
45 static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
46 static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
47 #else
48 static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
49 static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
50 static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
51 static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
52 #endif
53 #define pCf(z) (*_pCf(z))
54 #define pCd(z) (*_pCd(z))
55 typedef int logical;
56 typedef short int shortlogical;
57 typedef char logical1;
58 typedef char integer1;
59
60 #define TRUE_ (1)
61 #define FALSE_ (0)
62
63 /* Extern is for use with -E */
64 #ifndef Extern
65 #define Extern extern
66 #endif
67
68 /* I/O stuff */
69
70 typedef int flag;
71 typedef int ftnlen;
72 typedef int ftnint;
73
74 /*external read, write*/
75 typedef struct
76 {       flag cierr;
77         ftnint ciunit;
78         flag ciend;
79         char *cifmt;
80         ftnint cirec;
81 } cilist;
82
83 /*internal read, write*/
84 typedef struct
85 {       flag icierr;
86         char *iciunit;
87         flag iciend;
88         char *icifmt;
89         ftnint icirlen;
90         ftnint icirnum;
91 } icilist;
92
93 /*open*/
94 typedef struct
95 {       flag oerr;
96         ftnint ounit;
97         char *ofnm;
98         ftnlen ofnmlen;
99         char *osta;
100         char *oacc;
101         char *ofm;
102         ftnint orl;
103         char *oblnk;
104 } olist;
105
106 /*close*/
107 typedef struct
108 {       flag cerr;
109         ftnint cunit;
110         char *csta;
111 } cllist;
112
113 /*rewind, backspace, endfile*/
114 typedef struct
115 {       flag aerr;
116         ftnint aunit;
117 } alist;
118
119 /* inquire */
120 typedef struct
121 {       flag inerr;
122         ftnint inunit;
123         char *infile;
124         ftnlen infilen;
125         ftnint  *inex;  /*parameters in standard's order*/
126         ftnint  *inopen;
127         ftnint  *innum;
128         ftnint  *innamed;
129         char    *inname;
130         ftnlen  innamlen;
131         char    *inacc;
132         ftnlen  inacclen;
133         char    *inseq;
134         ftnlen  inseqlen;
135         char    *indir;
136         ftnlen  indirlen;
137         char    *infmt;
138         ftnlen  infmtlen;
139         char    *inform;
140         ftnint  informlen;
141         char    *inunf;
142         ftnlen  inunflen;
143         ftnint  *inrecl;
144         ftnint  *innrec;
145         char    *inblank;
146         ftnlen  inblanklen;
147 } inlist;
148
149 #define VOID void
150
151 union Multitype {       /* for multiple entry points */
152         integer1 g;
153         shortint h;
154         integer i;
155         /* longint j; */
156         real r;
157         doublereal d;
158         complex c;
159         doublecomplex z;
160         };
161
162 typedef union Multitype Multitype;
163
164 struct Vardesc {        /* for Namelist */
165         char *name;
166         char *addr;
167         ftnlen *dims;
168         int  type;
169         };
170 typedef struct Vardesc Vardesc;
171
172 struct Namelist {
173         char *name;
174         Vardesc **vars;
175         int nvars;
176         };
177 typedef struct Namelist Namelist;
178
179 #define abs(x) ((x) >= 0 ? (x) : -(x))
180 #define dabs(x) (fabs(x))
181 #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
182 #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
183 #define dmin(a,b) (f2cmin(a,b))
184 #define dmax(a,b) (f2cmax(a,b))
185 #define bit_test(a,b)   ((a) >> (b) & 1)
186 #define bit_clear(a,b)  ((a) & ~((uinteger)1 << (b)))
187 #define bit_set(a,b)    ((a) |  ((uinteger)1 << (b)))
188
189 #define abort_() { sig_die("Fortran abort routine called", 1); }
190 #define c_abs(z) (cabsf(Cf(z)))
191 #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
192 #ifdef _MSC_VER
193 #define c_div(c, a, b) {Cf(c)._Val[0] = (Cf(a)._Val[0]/Cf(b)._Val[0]); Cf(c)._Val[1]=(Cf(a)._Val[1]/Cf(b)._Val[1]);}
194 #define z_div(c, a, b) {Cd(c)._Val[0] = (Cd(a)._Val[0]/Cd(b)._Val[0]); Cd(c)._Val[1]=(Cd(a)._Val[1]/df(b)._Val[1]);}
195 #else
196 #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
197 #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
198 #endif
199 #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
200 #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
201 #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
202 //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
203 #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
204 #define d_abs(x) (fabs(*(x)))
205 #define d_acos(x) (acos(*(x)))
206 #define d_asin(x) (asin(*(x)))
207 #define d_atan(x) (atan(*(x)))
208 #define d_atn2(x, y) (atan2(*(x),*(y)))
209 #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
210 #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
211 #define d_cos(x) (cos(*(x)))
212 #define d_cosh(x) (cosh(*(x)))
213 #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
214 #define d_exp(x) (exp(*(x)))
215 #define d_imag(z) (cimag(Cd(z)))
216 #define r_imag(z) (cimagf(Cf(z)))
217 #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
218 #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
219 #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
220 #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
221 #define d_log(x) (log(*(x)))
222 #define d_mod(x, y) (fmod(*(x), *(y)))
223 #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
224 #define d_nint(x) u_nint(*(x))
225 #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
226 #define d_sign(a,b) u_sign(*(a),*(b))
227 #define r_sign(a,b) u_sign(*(a),*(b))
228 #define d_sin(x) (sin(*(x)))
229 #define d_sinh(x) (sinh(*(x)))
230 #define d_sqrt(x) (sqrt(*(x)))
231 #define d_tan(x) (tan(*(x)))
232 #define d_tanh(x) (tanh(*(x)))
233 #define i_abs(x) abs(*(x))
234 #define i_dnnt(x) ((integer)u_nint(*(x)))
235 #define i_len(s, n) (n)
236 #define i_nint(x) ((integer)u_nint(*(x)))
237 #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
238 #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
239 #define pow_si(B,E) spow_ui(*(B),*(E))
240 #define pow_ri(B,E) spow_ui(*(B),*(E))
241 #define pow_di(B,E) dpow_ui(*(B),*(E))
242 #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
243 #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
244 #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
245 #define s_cat(lpp, rpp, rnp, np, llp) {         ftnlen i, nc, ll; char *f__rp, *lp;     ll = (llp); lp = (lpp);         for(i=0; i < (int)*(np); ++i) {                 nc = ll;                if((rnp)[i] < nc) nc = (rnp)[i];                ll -= nc;               f__rp = (rpp)[i];               while(--nc >= 0) *lp++ = *(f__rp)++;         }  while(--ll >= 0) *lp++ = ' '; }
246 #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
247 #define s_copy(A,B,C,D) { int __i,__m; for (__i=0, __m=f2cmin((C),(D)); __i<__m && (B)[__i] != 0; ++__i) (A)[__i] = (B)[__i]; }
248 #define sig_die(s, kill) { exit(1); }
249 #define s_stop(s, n) {exit(0);}
250 static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
251 #define z_abs(z) (cabs(Cd(z)))
252 #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
253 #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
254 #define myexit_() break;
255 #define mycycle() continue;
256 #define myceiling(w) {ceil(w)}
257 #define myhuge(w) {HUGE_VAL}
258 //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
259 #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
260
261 /* procedure parameter types for -A and -C++ */
262
263 #define F2C_proc_par_types 1
264 #ifdef __cplusplus
265 typedef logical (*L_fp)(...);
266 #else
267 typedef logical (*L_fp)();
268 #endif
269
270 static float spow_ui(float x, integer n) {
271         float pow=1.0; unsigned long int u;
272         if(n != 0) {
273                 if(n < 0) n = -n, x = 1/x;
274                 for(u = n; ; ) {
275                         if(u & 01) pow *= x;
276                         if(u >>= 1) x *= x;
277                         else break;
278                 }
279         }
280         return pow;
281 }
282 static double dpow_ui(double x, integer n) {
283         double pow=1.0; unsigned long int u;
284         if(n != 0) {
285                 if(n < 0) n = -n, x = 1/x;
286                 for(u = n; ; ) {
287                         if(u & 01) pow *= x;
288                         if(u >>= 1) x *= x;
289                         else break;
290                 }
291         }
292         return pow;
293 }
294 #ifdef _MSC_VER
295 static _Fcomplex cpow_ui(complex x, integer n) {
296         complex pow={1.0,0.0}; unsigned long int u;
297                 if(n != 0) {
298                 if(n < 0) n = -n, x.r = 1/x.r, x.i=1/x.i;
299                 for(u = n; ; ) {
300                         if(u & 01) pow.r *= x.r, pow.i *= x.i;
301                         if(u >>= 1) x.r *= x.r, x.i *= x.i;
302                         else break;
303                 }
304         }
305         _Fcomplex p={pow.r, pow.i};
306         return p;
307 }
308 #else
309 static _Complex float cpow_ui(_Complex float x, integer n) {
310         _Complex float pow=1.0; unsigned long int u;
311         if(n != 0) {
312                 if(n < 0) n = -n, x = 1/x;
313                 for(u = n; ; ) {
314                         if(u & 01) pow *= x;
315                         if(u >>= 1) x *= x;
316                         else break;
317                 }
318         }
319         return pow;
320 }
321 #endif
322 #ifdef _MSC_VER
323 static _Dcomplex zpow_ui(_Dcomplex x, integer n) {
324         _Dcomplex pow={1.0,0.0}; unsigned long int u;
325         if(n != 0) {
326                 if(n < 0) n = -n, x._Val[0] = 1/x._Val[0], x._Val[1] =1/x._Val[1];
327                 for(u = n; ; ) {
328                         if(u & 01) pow._Val[0] *= x._Val[0], pow._Val[1] *= x._Val[1];
329                         if(u >>= 1) x._Val[0] *= x._Val[0], x._Val[1] *= x._Val[1];
330                         else break;
331                 }
332         }
333         _Dcomplex p = {pow._Val[0], pow._Val[1]};
334         return p;
335 }
336 #else
337 static _Complex double zpow_ui(_Complex double x, integer n) {
338         _Complex double pow=1.0; unsigned long int u;
339         if(n != 0) {
340                 if(n < 0) n = -n, x = 1/x;
341                 for(u = n; ; ) {
342                         if(u & 01) pow *= x;
343                         if(u >>= 1) x *= x;
344                         else break;
345                 }
346         }
347         return pow;
348 }
349 #endif
350 static integer pow_ii(integer x, integer n) {
351         integer pow; unsigned long int u;
352         if (n <= 0) {
353                 if (n == 0 || x == 1) pow = 1;
354                 else if (x != -1) pow = x == 0 ? 1/x : 0;
355                 else n = -n;
356         }
357         if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
358                 u = n;
359                 for(pow = 1; ; ) {
360                         if(u & 01) pow *= x;
361                         if(u >>= 1) x *= x;
362                         else break;
363                 }
364         }
365         return pow;
366 }
367 static integer dmaxloc_(double *w, integer s, integer e, integer *n)
368 {
369         double m; integer i, mi;
370         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
371                 if (w[i-1]>m) mi=i ,m=w[i-1];
372         return mi-s+1;
373 }
374 static integer smaxloc_(float *w, integer s, integer e, integer *n)
375 {
376         float m; integer i, mi;
377         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
378                 if (w[i-1]>m) mi=i ,m=w[i-1];
379         return mi-s+1;
380 }
381 static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
382         integer n = *n_, incx = *incx_, incy = *incy_, i;
383 #ifdef _MSC_VER
384         _Fcomplex zdotc = {0.0, 0.0};
385         if (incx == 1 && incy == 1) {
386                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
387                         zdotc._Val[0] += conjf(Cf(&x[i]))._Val[0] * Cf(&y[i])._Val[0];
388                         zdotc._Val[1] += conjf(Cf(&x[i]))._Val[1] * Cf(&y[i])._Val[1];
389                 }
390         } else {
391                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
392                         zdotc._Val[0] += conjf(Cf(&x[i*incx]))._Val[0] * Cf(&y[i*incy])._Val[0];
393                         zdotc._Val[1] += conjf(Cf(&x[i*incx]))._Val[1] * Cf(&y[i*incy])._Val[1];
394                 }
395         }
396         pCf(z) = zdotc;
397 }
398 #else
399         _Complex float zdotc = 0.0;
400         if (incx == 1 && incy == 1) {
401                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
402                         zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
403                 }
404         } else {
405                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
406                         zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
407                 }
408         }
409         pCf(z) = zdotc;
410 }
411 #endif
412 static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
413         integer n = *n_, incx = *incx_, incy = *incy_, i;
414 #ifdef _MSC_VER
415         _Dcomplex zdotc = {0.0, 0.0};
416         if (incx == 1 && incy == 1) {
417                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
418                         zdotc._Val[0] += conj(Cd(&x[i]))._Val[0] * Cd(&y[i])._Val[0];
419                         zdotc._Val[1] += conj(Cd(&x[i]))._Val[1] * Cd(&y[i])._Val[1];
420                 }
421         } else {
422                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
423                         zdotc._Val[0] += conj(Cd(&x[i*incx]))._Val[0] * Cd(&y[i*incy])._Val[0];
424                         zdotc._Val[1] += conj(Cd(&x[i*incx]))._Val[1] * Cd(&y[i*incy])._Val[1];
425                 }
426         }
427         pCd(z) = zdotc;
428 }
429 #else
430         _Complex double zdotc = 0.0;
431         if (incx == 1 && incy == 1) {
432                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
433                         zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
434                 }
435         } else {
436                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
437                         zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
438                 }
439         }
440         pCd(z) = zdotc;
441 }
442 #endif  
443 static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
444         integer n = *n_, incx = *incx_, incy = *incy_, i;
445 #ifdef _MSC_VER
446         _Fcomplex zdotc = {0.0, 0.0};
447         if (incx == 1 && incy == 1) {
448                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
449                         zdotc._Val[0] += Cf(&x[i])._Val[0] * Cf(&y[i])._Val[0];
450                         zdotc._Val[1] += Cf(&x[i])._Val[1] * Cf(&y[i])._Val[1];
451                 }
452         } else {
453                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
454                         zdotc._Val[0] += Cf(&x[i*incx])._Val[0] * Cf(&y[i*incy])._Val[0];
455                         zdotc._Val[1] += Cf(&x[i*incx])._Val[1] * Cf(&y[i*incy])._Val[1];
456                 }
457         }
458         pCf(z) = zdotc;
459 }
460 #else
461         _Complex float zdotc = 0.0;
462         if (incx == 1 && incy == 1) {
463                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
464                         zdotc += Cf(&x[i]) * Cf(&y[i]);
465                 }
466         } else {
467                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
468                         zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
469                 }
470         }
471         pCf(z) = zdotc;
472 }
473 #endif
474 static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
475         integer n = *n_, incx = *incx_, incy = *incy_, i;
476 #ifdef _MSC_VER
477         _Dcomplex zdotc = {0.0, 0.0};
478         if (incx == 1 && incy == 1) {
479                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
480                         zdotc._Val[0] += Cd(&x[i])._Val[0] * Cd(&y[i])._Val[0];
481                         zdotc._Val[1] += Cd(&x[i])._Val[1] * Cd(&y[i])._Val[1];
482                 }
483         } else {
484                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
485                         zdotc._Val[0] += Cd(&x[i*incx])._Val[0] * Cd(&y[i*incy])._Val[0];
486                         zdotc._Val[1] += Cd(&x[i*incx])._Val[1] * Cd(&y[i*incy])._Val[1];
487                 }
488         }
489         pCd(z) = zdotc;
490 }
491 #else
492         _Complex double zdotc = 0.0;
493         if (incx == 1 && incy == 1) {
494                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
495                         zdotc += Cd(&x[i]) * Cd(&y[i]);
496                 }
497         } else {
498                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
499                         zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
500                 }
501         }
502         pCd(z) = zdotc;
503 }
504 #endif
505 /*  -- translated by f2c (version 20000121).
506    You must link the resulting object file with the libraries:
507         -lf2c -lm   (in that order)
508 */
509
510
511
512
513 /* Table of constant values */
514
515 static integer c__1 = 1;
516 static integer c__0 = 0;
517 static integer c_n1 = -1;
518
519 /* > \brief <b> ZGEEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matr
520 ices</b> */
521
522 /*  =========== DOCUMENTATION =========== */
523
524 /* Online html documentation available at */
525 /*            http://www.netlib.org/lapack/explore-html/ */
526
527 /* > \htmlonly */
528 /* > Download ZGEEV + dependencies */
529 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgeev.f
530 "> */
531 /* > [TGZ]</a> */
532 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgeev.f
533 "> */
534 /* > [ZIP]</a> */
535 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgeev.f
536 "> */
537 /* > [TXT]</a> */
538 /* > \endhtmlonly */
539
540 /*  Definition: */
541 /*  =========== */
542
543 /*       SUBROUTINE ZGEEV( JOBVL, JOBVR, N, A, LDA, W, VL, LDVL, VR, LDVR, */
544 /*                         WORK, LWORK, RWORK, INFO ) */
545
546 /*       CHARACTER          JOBVL, JOBVR */
547 /*       INTEGER            INFO, LDA, LDVL, LDVR, LWORK, N */
548 /*       DOUBLE PRECISION   RWORK( * ) */
549 /*       COMPLEX*16         A( LDA, * ), VL( LDVL, * ), VR( LDVR, * ), */
550 /*      $                   W( * ), WORK( * ) */
551
552
553 /* > \par Purpose: */
554 /*  ============= */
555 /* > */
556 /* > \verbatim */
557 /* > */
558 /* > ZGEEV computes for an N-by-N complex nonsymmetric matrix A, the */
559 /* > eigenvalues and, optionally, the left and/or right eigenvectors. */
560 /* > */
561 /* > The right eigenvector v(j) of A satisfies */
562 /* >                  A * v(j) = lambda(j) * v(j) */
563 /* > where lambda(j) is its eigenvalue. */
564 /* > The left eigenvector u(j) of A satisfies */
565 /* >               u(j)**H * A = lambda(j) * u(j)**H */
566 /* > where u(j)**H denotes the conjugate transpose of u(j). */
567 /* > */
568 /* > The computed eigenvectors are normalized to have Euclidean norm */
569 /* > equal to 1 and largest component real. */
570 /* > \endverbatim */
571
572 /*  Arguments: */
573 /*  ========== */
574
575 /* > \param[in] JOBVL */
576 /* > \verbatim */
577 /* >          JOBVL is CHARACTER*1 */
578 /* >          = 'N': left eigenvectors of A are not computed; */
579 /* >          = 'V': left eigenvectors of are computed. */
580 /* > \endverbatim */
581 /* > */
582 /* > \param[in] JOBVR */
583 /* > \verbatim */
584 /* >          JOBVR is CHARACTER*1 */
585 /* >          = 'N': right eigenvectors of A are not computed; */
586 /* >          = 'V': right eigenvectors of A are computed. */
587 /* > \endverbatim */
588 /* > */
589 /* > \param[in] N */
590 /* > \verbatim */
591 /* >          N is INTEGER */
592 /* >          The order of the matrix A. N >= 0. */
593 /* > \endverbatim */
594 /* > */
595 /* > \param[in,out] A */
596 /* > \verbatim */
597 /* >          A is COMPLEX*16 array, dimension (LDA,N) */
598 /* >          On entry, the N-by-N matrix A. */
599 /* >          On exit, A has been overwritten. */
600 /* > \endverbatim */
601 /* > */
602 /* > \param[in] LDA */
603 /* > \verbatim */
604 /* >          LDA is INTEGER */
605 /* >          The leading dimension of the array A.  LDA >= f2cmax(1,N). */
606 /* > \endverbatim */
607 /* > */
608 /* > \param[out] W */
609 /* > \verbatim */
610 /* >          W is COMPLEX*16 array, dimension (N) */
611 /* >          W contains the computed eigenvalues. */
612 /* > \endverbatim */
613 /* > */
614 /* > \param[out] VL */
615 /* > \verbatim */
616 /* >          VL is COMPLEX*16 array, dimension (LDVL,N) */
617 /* >          If JOBVL = 'V', the left eigenvectors u(j) are stored one */
618 /* >          after another in the columns of VL, in the same order */
619 /* >          as their eigenvalues. */
620 /* >          If JOBVL = 'N', VL is not referenced. */
621 /* >          u(j) = VL(:,j), the j-th column of VL. */
622 /* > \endverbatim */
623 /* > */
624 /* > \param[in] LDVL */
625 /* > \verbatim */
626 /* >          LDVL is INTEGER */
627 /* >          The leading dimension of the array VL.  LDVL >= 1; if */
628 /* >          JOBVL = 'V', LDVL >= N. */
629 /* > \endverbatim */
630 /* > */
631 /* > \param[out] VR */
632 /* > \verbatim */
633 /* >          VR is COMPLEX*16 array, dimension (LDVR,N) */
634 /* >          If JOBVR = 'V', the right eigenvectors v(j) are stored one */
635 /* >          after another in the columns of VR, in the same order */
636 /* >          as their eigenvalues. */
637 /* >          If JOBVR = 'N', VR is not referenced. */
638 /* >          v(j) = VR(:,j), the j-th column of VR. */
639 /* > \endverbatim */
640 /* > */
641 /* > \param[in] LDVR */
642 /* > \verbatim */
643 /* >          LDVR is INTEGER */
644 /* >          The leading dimension of the array VR.  LDVR >= 1; if */
645 /* >          JOBVR = 'V', LDVR >= N. */
646 /* > \endverbatim */
647 /* > */
648 /* > \param[out] WORK */
649 /* > \verbatim */
650 /* >          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) */
651 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
652 /* > \endverbatim */
653 /* > */
654 /* > \param[in] LWORK */
655 /* > \verbatim */
656 /* >          LWORK is INTEGER */
657 /* >          The dimension of the array WORK.  LWORK >= f2cmax(1,2*N). */
658 /* >          For good performance, LWORK must generally be larger. */
659 /* > */
660 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
661 /* >          only calculates the optimal size of the WORK array, returns */
662 /* >          this value as the first entry of the WORK array, and no error */
663 /* >          message related to LWORK is issued by XERBLA. */
664 /* > \endverbatim */
665 /* > */
666 /* > \param[out] RWORK */
667 /* > \verbatim */
668 /* >          RWORK is DOUBLE PRECISION array, dimension (2*N) */
669 /* > \endverbatim */
670 /* > */
671 /* > \param[out] INFO */
672 /* > \verbatim */
673 /* >          INFO is INTEGER */
674 /* >          = 0:  successful exit */
675 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value. */
676 /* >          > 0:  if INFO = i, the QR algorithm failed to compute all the */
677 /* >                eigenvalues, and no eigenvectors have been computed; */
678 /* >                elements i+1:N of W contain eigenvalues which have */
679 /* >                converged. */
680 /* > \endverbatim */
681
682 /*  Authors: */
683 /*  ======== */
684
685 /* > \author Univ. of Tennessee */
686 /* > \author Univ. of California Berkeley */
687 /* > \author Univ. of Colorado Denver */
688 /* > \author NAG Ltd. */
689
690 /* > \date June 2016 */
691
692 /*  @precisions fortran z -> c */
693
694 /* > \ingroup complex16GEeigen */
695
696 /*  ===================================================================== */
697 /* Subroutine */ int zgeev_(char *jobvl, char *jobvr, integer *n, 
698         doublecomplex *a, integer *lda, doublecomplex *w, doublecomplex *vl, 
699         integer *ldvl, doublecomplex *vr, integer *ldvr, doublecomplex *work, 
700         integer *lwork, doublereal *rwork, integer *info)
701 {
702     /* System generated locals */
703     integer a_dim1, a_offset, vl_dim1, vl_offset, vr_dim1, vr_offset, i__1, 
704             i__2, i__3;
705     doublereal d__1, d__2;
706     doublecomplex z__1, z__2;
707
708     /* Local variables */
709     integer ibal;
710     char side[1];
711     doublereal anrm;
712     integer ierr, itau, iwrk, nout, i__, k;
713     extern logical lsame_(char *, char *);
714     extern /* Subroutine */ int zscal_(integer *, doublecomplex *, 
715             doublecomplex *, integer *), dlabad_(doublereal *, doublereal *);
716     extern doublereal dznrm2_(integer *, doublecomplex *, integer *);
717     logical scalea;
718     extern doublereal dlamch_(char *);
719     doublereal cscale;
720     extern /* Subroutine */ int zgebak_(char *, char *, integer *, integer *, 
721             integer *, doublereal *, integer *, doublecomplex *, integer *, 
722             integer *), zgebal_(char *, integer *, 
723             doublecomplex *, integer *, integer *, integer *, doublereal *, 
724             integer *);
725     extern integer idamax_(integer *, doublereal *, integer *);
726     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
727     extern integer ilaenv_(integer *, char *, char *, integer *, integer *, 
728             integer *, integer *, ftnlen, ftnlen);
729     logical select[1];
730     extern /* Subroutine */ int zdscal_(integer *, doublereal *, 
731             doublecomplex *, integer *);
732     doublereal bignum;
733     extern doublereal zlange_(char *, integer *, integer *, doublecomplex *, 
734             integer *, doublereal *);
735     extern /* Subroutine */ int zgehrd_(integer *, integer *, integer *, 
736             doublecomplex *, integer *, doublecomplex *, doublecomplex *, 
737             integer *, integer *), zlascl_(char *, integer *, integer *, 
738             doublereal *, doublereal *, integer *, integer *, doublecomplex *,
739              integer *, integer *), zlacpy_(char *, integer *, 
740             integer *, doublecomplex *, integer *, doublecomplex *, integer *);
741     integer minwrk, maxwrk;
742     logical wantvl;
743     doublereal smlnum;
744     integer hswork, irwork;
745     extern /* Subroutine */ int zhseqr_(char *, char *, integer *, integer *, 
746             integer *, doublecomplex *, integer *, doublecomplex *, 
747             doublecomplex *, integer *, doublecomplex *, integer *, integer *), zunghr_(integer *, integer *, integer *, 
748             doublecomplex *, integer *, doublecomplex *, doublecomplex *, 
749             integer *, integer *);
750     logical lquery, wantvr;
751     integer ihi;
752     extern /* Subroutine */ int ztrevc3_(char *, char *, logical *, integer *,
753              doublecomplex *, integer *, doublecomplex *, integer *, 
754             doublecomplex *, integer *, integer *, integer *, doublecomplex *,
755              integer *, doublereal *, integer *, integer *);
756     doublereal scl;
757     integer ilo;
758     doublereal dum[1], eps;
759     doublecomplex tmp;
760     integer lwork_trevc__;
761
762
763 /*  -- LAPACK driver routine (version 3.7.0) -- */
764 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
765 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
766 /*     June 2016 */
767
768
769 /*  ===================================================================== */
770
771
772 /*     Test the input arguments */
773
774     /* Parameter adjustments */
775     a_dim1 = *lda;
776     a_offset = 1 + a_dim1 * 1;
777     a -= a_offset;
778     --w;
779     vl_dim1 = *ldvl;
780     vl_offset = 1 + vl_dim1 * 1;
781     vl -= vl_offset;
782     vr_dim1 = *ldvr;
783     vr_offset = 1 + vr_dim1 * 1;
784     vr -= vr_offset;
785     --work;
786     --rwork;
787
788     /* Function Body */
789     *info = 0;
790     lquery = *lwork == -1;
791     wantvl = lsame_(jobvl, "V");
792     wantvr = lsame_(jobvr, "V");
793     if (! wantvl && ! lsame_(jobvl, "N")) {
794         *info = -1;
795     } else if (! wantvr && ! lsame_(jobvr, "N")) {
796         *info = -2;
797     } else if (*n < 0) {
798         *info = -3;
799     } else if (*lda < f2cmax(1,*n)) {
800         *info = -5;
801     } else if (*ldvl < 1 || wantvl && *ldvl < *n) {
802         *info = -8;
803     } else if (*ldvr < 1 || wantvr && *ldvr < *n) {
804         *info = -10;
805     }
806
807 /*     Compute workspace */
808 /*      (Note: Comments in the code beginning "Workspace:" describe the */
809 /*       minimal amount of workspace needed at that point in the code, */
810 /*       as well as the preferred amount for good performance. */
811 /*       CWorkspace refers to complex workspace, and RWorkspace to real */
812 /*       workspace. NB refers to the optimal block size for the */
813 /*       immediately following subroutine, as returned by ILAENV. */
814 /*       HSWORK refers to the workspace preferred by ZHSEQR, as */
815 /*       calculated below. HSWORK is computed assuming ILO=1 and IHI=N, */
816 /*       the worst case.) */
817
818     if (*info == 0) {
819         if (*n == 0) {
820             minwrk = 1;
821             maxwrk = 1;
822         } else {
823             maxwrk = *n + *n * ilaenv_(&c__1, "ZGEHRD", " ", n, &c__1, n, &
824                     c__0, (ftnlen)6, (ftnlen)1);
825             minwrk = *n << 1;
826             if (wantvl) {
827 /* Computing MAX */
828                 i__1 = maxwrk, i__2 = *n + (*n - 1) * ilaenv_(&c__1, "ZUNGHR",
829                          " ", n, &c__1, n, &c_n1, (ftnlen)6, (ftnlen)1);
830                 maxwrk = f2cmax(i__1,i__2);
831                 ztrevc3_("L", "B", select, n, &a[a_offset], lda, &vl[
832                         vl_offset], ldvl, &vr[vr_offset], ldvr, n, &nout, &
833                         work[1], &c_n1, &rwork[1], &c_n1, &ierr);
834                 lwork_trevc__ = (integer) work[1].r;
835 /* Computing MAX */
836                 i__1 = maxwrk, i__2 = *n + lwork_trevc__;
837                 maxwrk = f2cmax(i__1,i__2);
838                 zhseqr_("S", "V", n, &c__1, n, &a[a_offset], lda, &w[1], &vl[
839                         vl_offset], ldvl, &work[1], &c_n1, info);
840             } else if (wantvr) {
841 /* Computing MAX */
842                 i__1 = maxwrk, i__2 = *n + (*n - 1) * ilaenv_(&c__1, "ZUNGHR",
843                          " ", n, &c__1, n, &c_n1, (ftnlen)6, (ftnlen)1);
844                 maxwrk = f2cmax(i__1,i__2);
845                 ztrevc3_("R", "B", select, n, &a[a_offset], lda, &vl[
846                         vl_offset], ldvl, &vr[vr_offset], ldvr, n, &nout, &
847                         work[1], &c_n1, &rwork[1], &c_n1, &ierr);
848                 lwork_trevc__ = (integer) work[1].r;
849 /* Computing MAX */
850                 i__1 = maxwrk, i__2 = *n + lwork_trevc__;
851                 maxwrk = f2cmax(i__1,i__2);
852                 zhseqr_("S", "V", n, &c__1, n, &a[a_offset], lda, &w[1], &vr[
853                         vr_offset], ldvr, &work[1], &c_n1, info);
854             } else {
855                 zhseqr_("E", "N", n, &c__1, n, &a[a_offset], lda, &w[1], &vr[
856                         vr_offset], ldvr, &work[1], &c_n1, info);
857             }
858             hswork = (integer) work[1].r;
859 /* Computing MAX */
860             i__1 = f2cmax(maxwrk,hswork);
861             maxwrk = f2cmax(i__1,minwrk);
862         }
863         work[1].r = (doublereal) maxwrk, work[1].i = 0.;
864
865         if (*lwork < minwrk && ! lquery) {
866             *info = -12;
867         }
868     }
869
870     if (*info != 0) {
871         i__1 = -(*info);
872         xerbla_("ZGEEV ", &i__1, (ftnlen)6);
873         return 0;
874     } else if (lquery) {
875         return 0;
876     }
877
878 /*     Quick return if possible */
879
880     if (*n == 0) {
881         return 0;
882     }
883
884 /*     Get machine constants */
885
886     eps = dlamch_("P");
887     smlnum = dlamch_("S");
888     bignum = 1. / smlnum;
889     dlabad_(&smlnum, &bignum);
890     smlnum = sqrt(smlnum) / eps;
891     bignum = 1. / smlnum;
892
893 /*     Scale A if f2cmax element outside range [SMLNUM,BIGNUM] */
894
895     anrm = zlange_("M", n, n, &a[a_offset], lda, dum);
896     scalea = FALSE_;
897     if (anrm > 0. && anrm < smlnum) {
898         scalea = TRUE_;
899         cscale = smlnum;
900     } else if (anrm > bignum) {
901         scalea = TRUE_;
902         cscale = bignum;
903     }
904     if (scalea) {
905         zlascl_("G", &c__0, &c__0, &anrm, &cscale, n, n, &a[a_offset], lda, &
906                 ierr);
907     }
908
909 /*     Balance the matrix */
910 /*     (CWorkspace: none) */
911 /*     (RWorkspace: need N) */
912
913     ibal = 1;
914     zgebal_("B", n, &a[a_offset], lda, &ilo, &ihi, &rwork[ibal], &ierr);
915
916 /*     Reduce to upper Hessenberg form */
917 /*     (CWorkspace: need 2*N, prefer N+N*NB) */
918 /*     (RWorkspace: none) */
919
920     itau = 1;
921     iwrk = itau + *n;
922     i__1 = *lwork - iwrk + 1;
923     zgehrd_(n, &ilo, &ihi, &a[a_offset], lda, &work[itau], &work[iwrk], &i__1,
924              &ierr);
925
926     if (wantvl) {
927
928 /*        Want left eigenvectors */
929 /*        Copy Householder vectors to VL */
930
931         *(unsigned char *)side = 'L';
932         zlacpy_("L", n, n, &a[a_offset], lda, &vl[vl_offset], ldvl)
933                 ;
934
935 /*        Generate unitary matrix in VL */
936 /*        (CWorkspace: need 2*N-1, prefer N+(N-1)*NB) */
937 /*        (RWorkspace: none) */
938
939         i__1 = *lwork - iwrk + 1;
940         zunghr_(n, &ilo, &ihi, &vl[vl_offset], ldvl, &work[itau], &work[iwrk],
941                  &i__1, &ierr);
942
943 /*        Perform QR iteration, accumulating Schur vectors in VL */
944 /*        (CWorkspace: need 1, prefer HSWORK (see comments) ) */
945 /*        (RWorkspace: none) */
946
947         iwrk = itau;
948         i__1 = *lwork - iwrk + 1;
949         zhseqr_("S", "V", n, &ilo, &ihi, &a[a_offset], lda, &w[1], &vl[
950                 vl_offset], ldvl, &work[iwrk], &i__1, info);
951
952         if (wantvr) {
953
954 /*           Want left and right eigenvectors */
955 /*           Copy Schur vectors to VR */
956
957             *(unsigned char *)side = 'B';
958             zlacpy_("F", n, n, &vl[vl_offset], ldvl, &vr[vr_offset], ldvr);
959         }
960
961     } else if (wantvr) {
962
963 /*        Want right eigenvectors */
964 /*        Copy Householder vectors to VR */
965
966         *(unsigned char *)side = 'R';
967         zlacpy_("L", n, n, &a[a_offset], lda, &vr[vr_offset], ldvr)
968                 ;
969
970 /*        Generate unitary matrix in VR */
971 /*        (CWorkspace: need 2*N-1, prefer N+(N-1)*NB) */
972 /*        (RWorkspace: none) */
973
974         i__1 = *lwork - iwrk + 1;
975         zunghr_(n, &ilo, &ihi, &vr[vr_offset], ldvr, &work[itau], &work[iwrk],
976                  &i__1, &ierr);
977
978 /*        Perform QR iteration, accumulating Schur vectors in VR */
979 /*        (CWorkspace: need 1, prefer HSWORK (see comments) ) */
980 /*        (RWorkspace: none) */
981
982         iwrk = itau;
983         i__1 = *lwork - iwrk + 1;
984         zhseqr_("S", "V", n, &ilo, &ihi, &a[a_offset], lda, &w[1], &vr[
985                 vr_offset], ldvr, &work[iwrk], &i__1, info);
986
987     } else {
988
989 /*        Compute eigenvalues only */
990 /*        (CWorkspace: need 1, prefer HSWORK (see comments) ) */
991 /*        (RWorkspace: none) */
992
993         iwrk = itau;
994         i__1 = *lwork - iwrk + 1;
995         zhseqr_("E", "N", n, &ilo, &ihi, &a[a_offset], lda, &w[1], &vr[
996                 vr_offset], ldvr, &work[iwrk], &i__1, info);
997     }
998
999 /*     If INFO .NE. 0 from ZHSEQR, then quit */
1000
1001     if (*info != 0) {
1002         goto L50;
1003     }
1004
1005     if (wantvl || wantvr) {
1006
1007 /*        Compute left and/or right eigenvectors */
1008 /*        (CWorkspace: need 2*N, prefer N + 2*N*NB) */
1009 /*        (RWorkspace: need 2*N) */
1010
1011         irwork = ibal + *n;
1012         i__1 = *lwork - iwrk + 1;
1013         ztrevc3_(side, "B", select, n, &a[a_offset], lda, &vl[vl_offset], 
1014                 ldvl, &vr[vr_offset], ldvr, n, &nout, &work[iwrk], &i__1, &
1015                 rwork[irwork], n, &ierr);
1016     }
1017
1018     if (wantvl) {
1019
1020 /*        Undo balancing of left eigenvectors */
1021 /*        (CWorkspace: none) */
1022 /*        (RWorkspace: need N) */
1023
1024         zgebak_("B", "L", n, &ilo, &ihi, &rwork[ibal], n, &vl[vl_offset], 
1025                 ldvl, &ierr);
1026
1027 /*        Normalize left eigenvectors and make largest component real */
1028
1029         i__1 = *n;
1030         for (i__ = 1; i__ <= i__1; ++i__) {
1031             scl = 1. / dznrm2_(n, &vl[i__ * vl_dim1 + 1], &c__1);
1032             zdscal_(n, &scl, &vl[i__ * vl_dim1 + 1], &c__1);
1033             i__2 = *n;
1034             for (k = 1; k <= i__2; ++k) {
1035                 i__3 = k + i__ * vl_dim1;
1036 /* Computing 2nd power */
1037                 d__1 = vl[i__3].r;
1038 /* Computing 2nd power */
1039                 d__2 = d_imag(&vl[k + i__ * vl_dim1]);
1040                 rwork[irwork + k - 1] = d__1 * d__1 + d__2 * d__2;
1041 /*     $                               AIMAG( VL( K, I ) )**2 */
1042 /* L10: */
1043             }
1044             k = idamax_(n, &rwork[irwork], &c__1);
1045             d_cnjg(&z__2, &vl[k + i__ * vl_dim1]);
1046             d__1 = sqrt(rwork[irwork + k - 1]);
1047             z__1.r = z__2.r / d__1, z__1.i = z__2.i / d__1;
1048             tmp.r = z__1.r, tmp.i = z__1.i;
1049             zscal_(n, &tmp, &vl[i__ * vl_dim1 + 1], &c__1);
1050             i__2 = k + i__ * vl_dim1;
1051             i__3 = k + i__ * vl_dim1;
1052             d__1 = vl[i__3].r;
1053             z__1.r = d__1, z__1.i = 0.;
1054             vl[i__2].r = z__1.r, vl[i__2].i = z__1.i;
1055 /* L20: */
1056         }
1057     }
1058
1059     if (wantvr) {
1060
1061 /*        Undo balancing of right eigenvectors */
1062 /*        (CWorkspace: none) */
1063 /*        (RWorkspace: need N) */
1064
1065         zgebak_("B", "R", n, &ilo, &ihi, &rwork[ibal], n, &vr[vr_offset], 
1066                 ldvr, &ierr);
1067
1068 /*        Normalize right eigenvectors and make largest component real */
1069
1070         i__1 = *n;
1071         for (i__ = 1; i__ <= i__1; ++i__) {
1072             scl = 1. / dznrm2_(n, &vr[i__ * vr_dim1 + 1], &c__1);
1073             zdscal_(n, &scl, &vr[i__ * vr_dim1 + 1], &c__1);
1074             i__2 = *n;
1075             for (k = 1; k <= i__2; ++k) {
1076                 i__3 = k + i__ * vr_dim1;
1077 /* Computing 2nd power */
1078                 d__1 = vr[i__3].r;
1079 /* Computing 2nd power */
1080                 d__2 = d_imag(&vr[k + i__ * vr_dim1]);
1081                 rwork[irwork + k - 1] = d__1 * d__1 + d__2 * d__2;
1082 /*     $                               AIMAG( VR( K, I ) )**2 */
1083 /* L30: */
1084             }
1085             k = idamax_(n, &rwork[irwork], &c__1);
1086             d_cnjg(&z__2, &vr[k + i__ * vr_dim1]);
1087             d__1 = sqrt(rwork[irwork + k - 1]);
1088             z__1.r = z__2.r / d__1, z__1.i = z__2.i / d__1;
1089             tmp.r = z__1.r, tmp.i = z__1.i;
1090             zscal_(n, &tmp, &vr[i__ * vr_dim1 + 1], &c__1);
1091             i__2 = k + i__ * vr_dim1;
1092             i__3 = k + i__ * vr_dim1;
1093             d__1 = vr[i__3].r;
1094             z__1.r = d__1, z__1.i = 0.;
1095             vr[i__2].r = z__1.r, vr[i__2].i = z__1.i;
1096 /* L40: */
1097         }
1098     }
1099
1100 /*     Undo scaling if necessary */
1101
1102 L50:
1103     if (scalea) {
1104         i__1 = *n - *info;
1105 /* Computing MAX */
1106         i__3 = *n - *info;
1107         i__2 = f2cmax(i__3,1);
1108         zlascl_("G", &c__0, &c__0, &cscale, &anrm, &i__1, &c__1, &w[*info + 1]
1109                 , &i__2, &ierr);
1110         if (*info > 0) {
1111             i__1 = ilo - 1;
1112             zlascl_("G", &c__0, &c__0, &cscale, &anrm, &i__1, &c__1, &w[1], n,
1113                      &ierr);
1114         }
1115     }
1116
1117     work[1].r = (doublereal) maxwrk, work[1].i = 0.;
1118     return 0;
1119
1120 /*     End of ZGEEV */
1121
1122 } /* zgeev_ */
1123