C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / cgeevx.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
514 /* Table of constant values */
515
516 static integer c__1 = 1;
517 static integer c__0 = 0;
518 static integer c_n1 = -1;
519
520 /* > \brief <b> CGEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE mat
521 rices</b> */
522
523 /*  =========== DOCUMENTATION =========== */
524
525 /* Online html documentation available at */
526 /*            http://www.netlib.org/lapack/explore-html/ */
527
528 /* > \htmlonly */
529 /* > Download CGEEVX + dependencies */
530 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cgeevx.
531 f"> */
532 /* > [TGZ]</a> */
533 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cgeevx.
534 f"> */
535 /* > [ZIP]</a> */
536 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cgeevx.
537 f"> */
538 /* > [TXT]</a> */
539 /* > \endhtmlonly */
540
541 /*  Definition: */
542 /*  =========== */
543
544 /*       SUBROUTINE CGEEVX( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, W, VL, */
545 /*                          LDVL, VR, LDVR, ILO, IHI, SCALE, ABNRM, RCONDE, */
546 /*                          RCONDV, WORK, LWORK, RWORK, INFO ) */
547
548 /*       CHARACTER          BALANC, JOBVL, JOBVR, SENSE */
549 /*       INTEGER            IHI, ILO, INFO, LDA, LDVL, LDVR, LWORK, N */
550 /*       REAL               ABNRM */
551 /*       REAL               RCONDE( * ), RCONDV( * ), RWORK( * ), */
552 /*      $                   SCALE( * ) */
553 /*       COMPLEX            A( LDA, * ), VL( LDVL, * ), VR( LDVR, * ), */
554 /*      $                   W( * ), WORK( * ) */
555
556
557 /* > \par Purpose: */
558 /*  ============= */
559 /* > */
560 /* > \verbatim */
561 /* > */
562 /* > CGEEVX computes for an N-by-N complex nonsymmetric matrix A, the */
563 /* > eigenvalues and, optionally, the left and/or right eigenvectors. */
564 /* > */
565 /* > Optionally also, it computes a balancing transformation to improve */
566 /* > the conditioning of the eigenvalues and eigenvectors (ILO, IHI, */
567 /* > SCALE, and ABNRM), reciprocal condition numbers for the eigenvalues */
568 /* > (RCONDE), and reciprocal condition numbers for the right */
569 /* > eigenvectors (RCONDV). */
570 /* > */
571 /* > The right eigenvector v(j) of A satisfies */
572 /* >                  A * v(j) = lambda(j) * v(j) */
573 /* > where lambda(j) is its eigenvalue. */
574 /* > The left eigenvector u(j) of A satisfies */
575 /* >               u(j)**H * A = lambda(j) * u(j)**H */
576 /* > where u(j)**H denotes the conjugate transpose of u(j). */
577 /* > */
578 /* > The computed eigenvectors are normalized to have Euclidean norm */
579 /* > equal to 1 and largest component real. */
580 /* > */
581 /* > Balancing a matrix means permuting the rows and columns to make it */
582 /* > more nearly upper triangular, and applying a diagonal similarity */
583 /* > transformation D * A * D**(-1), where D is a diagonal matrix, to */
584 /* > make its rows and columns closer in norm and the condition numbers */
585 /* > of its eigenvalues and eigenvectors smaller.  The computed */
586 /* > reciprocal condition numbers correspond to the balanced matrix. */
587 /* > Permuting rows and columns will not change the condition numbers */
588 /* > (in exact arithmetic) but diagonal scaling will.  For further */
589 /* > explanation of balancing, see section 4.10.2 of the LAPACK */
590 /* > Users' Guide. */
591 /* > \endverbatim */
592
593 /*  Arguments: */
594 /*  ========== */
595
596 /* > \param[in] BALANC */
597 /* > \verbatim */
598 /* >          BALANC is CHARACTER*1 */
599 /* >          Indicates how the input matrix should be diagonally scaled */
600 /* >          and/or permuted to improve the conditioning of its */
601 /* >          eigenvalues. */
602 /* >          = 'N': Do not diagonally scale or permute; */
603 /* >          = 'P': Perform permutations to make the matrix more nearly */
604 /* >                 upper triangular. Do not diagonally scale; */
605 /* >          = 'S': Diagonally scale the matrix, ie. replace A by */
606 /* >                 D*A*D**(-1), where D is a diagonal matrix chosen */
607 /* >                 to make the rows and columns of A more equal in */
608 /* >                 norm. Do not permute; */
609 /* >          = 'B': Both diagonally scale and permute A. */
610 /* > */
611 /* >          Computed reciprocal condition numbers will be for the matrix */
612 /* >          after balancing and/or permuting. Permuting does not change */
613 /* >          condition numbers (in exact arithmetic), but balancing does. */
614 /* > \endverbatim */
615 /* > */
616 /* > \param[in] JOBVL */
617 /* > \verbatim */
618 /* >          JOBVL is CHARACTER*1 */
619 /* >          = 'N': left eigenvectors of A are not computed; */
620 /* >          = 'V': left eigenvectors of A are computed. */
621 /* >          If SENSE = 'E' or 'B', JOBVL must = 'V'. */
622 /* > \endverbatim */
623 /* > */
624 /* > \param[in] JOBVR */
625 /* > \verbatim */
626 /* >          JOBVR is CHARACTER*1 */
627 /* >          = 'N': right eigenvectors of A are not computed; */
628 /* >          = 'V': right eigenvectors of A are computed. */
629 /* >          If SENSE = 'E' or 'B', JOBVR must = 'V'. */
630 /* > \endverbatim */
631 /* > */
632 /* > \param[in] SENSE */
633 /* > \verbatim */
634 /* >          SENSE is CHARACTER*1 */
635 /* >          Determines which reciprocal condition numbers are computed. */
636 /* >          = 'N': None are computed; */
637 /* >          = 'E': Computed for eigenvalues only; */
638 /* >          = 'V': Computed for right eigenvectors only; */
639 /* >          = 'B': Computed for eigenvalues and right eigenvectors. */
640 /* > */
641 /* >          If SENSE = 'E' or 'B', both left and right eigenvectors */
642 /* >          must also be computed (JOBVL = 'V' and JOBVR = 'V'). */
643 /* > \endverbatim */
644 /* > */
645 /* > \param[in] N */
646 /* > \verbatim */
647 /* >          N is INTEGER */
648 /* >          The order of the matrix A. N >= 0. */
649 /* > \endverbatim */
650 /* > */
651 /* > \param[in,out] A */
652 /* > \verbatim */
653 /* >          A is COMPLEX array, dimension (LDA,N) */
654 /* >          On entry, the N-by-N matrix A. */
655 /* >          On exit, A has been overwritten.  If JOBVL = 'V' or */
656 /* >          JOBVR = 'V', A contains the Schur form of the balanced */
657 /* >          version of the matrix A. */
658 /* > \endverbatim */
659 /* > */
660 /* > \param[in] LDA */
661 /* > \verbatim */
662 /* >          LDA is INTEGER */
663 /* >          The leading dimension of the array A.  LDA >= f2cmax(1,N). */
664 /* > \endverbatim */
665 /* > */
666 /* > \param[out] W */
667 /* > \verbatim */
668 /* >          W is COMPLEX array, dimension (N) */
669 /* >          W contains the computed eigenvalues. */
670 /* > \endverbatim */
671 /* > */
672 /* > \param[out] VL */
673 /* > \verbatim */
674 /* >          VL is COMPLEX array, dimension (LDVL,N) */
675 /* >          If JOBVL = 'V', the left eigenvectors u(j) are stored one */
676 /* >          after another in the columns of VL, in the same order */
677 /* >          as their eigenvalues. */
678 /* >          If JOBVL = 'N', VL is not referenced. */
679 /* >          u(j) = VL(:,j), the j-th column of VL. */
680 /* > \endverbatim */
681 /* > */
682 /* > \param[in] LDVL */
683 /* > \verbatim */
684 /* >          LDVL is INTEGER */
685 /* >          The leading dimension of the array VL.  LDVL >= 1; if */
686 /* >          JOBVL = 'V', LDVL >= N. */
687 /* > \endverbatim */
688 /* > */
689 /* > \param[out] VR */
690 /* > \verbatim */
691 /* >          VR is COMPLEX array, dimension (LDVR,N) */
692 /* >          If JOBVR = 'V', the right eigenvectors v(j) are stored one */
693 /* >          after another in the columns of VR, in the same order */
694 /* >          as their eigenvalues. */
695 /* >          If JOBVR = 'N', VR is not referenced. */
696 /* >          v(j) = VR(:,j), the j-th column of VR. */
697 /* > \endverbatim */
698 /* > */
699 /* > \param[in] LDVR */
700 /* > \verbatim */
701 /* >          LDVR is INTEGER */
702 /* >          The leading dimension of the array VR.  LDVR >= 1; if */
703 /* >          JOBVR = 'V', LDVR >= N. */
704 /* > \endverbatim */
705 /* > */
706 /* > \param[out] ILO */
707 /* > \verbatim */
708 /* >          ILO is INTEGER */
709 /* > \endverbatim */
710 /* > */
711 /* > \param[out] IHI */
712 /* > \verbatim */
713 /* >          IHI is INTEGER */
714 /* >          ILO and IHI are integer values determined when A was */
715 /* >          balanced.  The balanced A(i,j) = 0 if I > J and */
716 /* >          J = 1,...,ILO-1 or I = IHI+1,...,N. */
717 /* > \endverbatim */
718 /* > */
719 /* > \param[out] SCALE */
720 /* > \verbatim */
721 /* >          SCALE is REAL array, dimension (N) */
722 /* >          Details of the permutations and scaling factors applied */
723 /* >          when balancing A.  If P(j) is the index of the row and column */
724 /* >          interchanged with row and column j, and D(j) is the scaling */
725 /* >          factor applied to row and column j, then */
726 /* >          SCALE(J) = P(J),    for J = 1,...,ILO-1 */
727 /* >                   = D(J),    for J = ILO,...,IHI */
728 /* >                   = P(J)     for J = IHI+1,...,N. */
729 /* >          The order in which the interchanges are made is N to IHI+1, */
730 /* >          then 1 to ILO-1. */
731 /* > \endverbatim */
732 /* > */
733 /* > \param[out] ABNRM */
734 /* > \verbatim */
735 /* >          ABNRM is REAL */
736 /* >          The one-norm of the balanced matrix (the maximum */
737 /* >          of the sum of absolute values of elements of any column). */
738 /* > \endverbatim */
739 /* > */
740 /* > \param[out] RCONDE */
741 /* > \verbatim */
742 /* >          RCONDE is REAL array, dimension (N) */
743 /* >          RCONDE(j) is the reciprocal condition number of the j-th */
744 /* >          eigenvalue. */
745 /* > \endverbatim */
746 /* > */
747 /* > \param[out] RCONDV */
748 /* > \verbatim */
749 /* >          RCONDV is REAL array, dimension (N) */
750 /* >          RCONDV(j) is the reciprocal condition number of the j-th */
751 /* >          right eigenvector. */
752 /* > \endverbatim */
753 /* > */
754 /* > \param[out] WORK */
755 /* > \verbatim */
756 /* >          WORK is COMPLEX array, dimension (MAX(1,LWORK)) */
757 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
758 /* > \endverbatim */
759 /* > */
760 /* > \param[in] LWORK */
761 /* > \verbatim */
762 /* >          LWORK is INTEGER */
763 /* >          The dimension of the array WORK.  If SENSE = 'N' or 'E', */
764 /* >          LWORK >= f2cmax(1,2*N), and if SENSE = 'V' or 'B', */
765 /* >          LWORK >= N*N+2*N. */
766 /* >          For good performance, LWORK must generally be larger. */
767 /* > */
768 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
769 /* >          only calculates the optimal size of the WORK array, returns */
770 /* >          this value as the first entry of the WORK array, and no error */
771 /* >          message related to LWORK is issued by XERBLA. */
772 /* > \endverbatim */
773 /* > */
774 /* > \param[out] RWORK */
775 /* > \verbatim */
776 /* >          RWORK is REAL array, dimension (2*N) */
777 /* > \endverbatim */
778 /* > */
779 /* > \param[out] INFO */
780 /* > \verbatim */
781 /* >          INFO is INTEGER */
782 /* >          = 0:  successful exit */
783 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value. */
784 /* >          > 0:  if INFO = i, the QR algorithm failed to compute all the */
785 /* >                eigenvalues, and no eigenvectors or condition numbers */
786 /* >                have been computed; elements 1:ILO-1 and i+1:N of W */
787 /* >                contain eigenvalues which have converged. */
788 /* > \endverbatim */
789
790 /*  Authors: */
791 /*  ======== */
792
793 /* > \author Univ. of Tennessee */
794 /* > \author Univ. of California Berkeley */
795 /* > \author Univ. of Colorado Denver */
796 /* > \author NAG Ltd. */
797
798 /* > \date June 2016 */
799
800 /*  @generated from zgeevx.f, fortran z -> c, Tue Apr 19 01:47:44 2016 */
801
802 /* > \ingroup complexGEeigen */
803
804 /*  ===================================================================== */
805 /* Subroutine */ int cgeevx_(char *balanc, char *jobvl, char *jobvr, char *
806         sense, integer *n, complex *a, integer *lda, complex *w, complex *vl, 
807         integer *ldvl, complex *vr, integer *ldvr, integer *ilo, integer *ihi,
808          real *scale, real *abnrm, real *rconde, real *rcondv, complex *work, 
809         integer *lwork, real *rwork, integer *info)
810 {
811     /* System generated locals */
812     integer a_dim1, a_offset, vl_dim1, vl_offset, vr_dim1, vr_offset, i__1, 
813             i__2, i__3;
814     real r__1, r__2;
815     complex q__1, q__2;
816
817     /* Local variables */
818     char side[1];
819     real anrm;
820     integer ierr, itau, iwrk, nout, i__, k;
821     extern /* Subroutine */ int cscal_(integer *, complex *, complex *, 
822             integer *);
823     integer icond;
824     extern logical lsame_(char *, char *);
825     extern real scnrm2_(integer *, complex *, integer *);
826     extern /* Subroutine */ int cgebak_(char *, char *, integer *, integer *, 
827             integer *, real *, integer *, complex *, integer *, integer *), cgebal_(char *, integer *, complex *, integer *, 
828             integer *, integer *, real *, integer *), slabad_(real *, 
829             real *);
830     logical scalea;
831     extern real clange_(char *, integer *, integer *, complex *, integer *, 
832             real *);
833     real cscale;
834     extern /* Subroutine */ int cgehrd_(integer *, integer *, integer *, 
835             complex *, integer *, complex *, complex *, integer *, integer *),
836              clascl_(char *, integer *, integer *, real *, real *, integer *, 
837             integer *, complex *, integer *, integer *);
838     extern real slamch_(char *);
839     extern /* Subroutine */ int csscal_(integer *, real *, complex *, integer 
840             *), clacpy_(char *, integer *, integer *, complex *, integer *, 
841             complex *, integer *), xerbla_(char *, integer *, ftnlen);
842     extern integer ilaenv_(integer *, char *, char *, integer *, integer *, 
843             integer *, integer *, ftnlen, ftnlen);
844     logical select[1];
845     real bignum;
846     extern /* Subroutine */ int slascl_(char *, integer *, integer *, real *, 
847             real *, integer *, integer *, real *, integer *, integer *);
848     extern integer isamax_(integer *, real *, integer *);
849     extern /* Subroutine */ int chseqr_(char *, char *, integer *, integer *, 
850             integer *, complex *, integer *, complex *, complex *, integer *, 
851             complex *, integer *, integer *), cunghr_(integer 
852             *, integer *, integer *, complex *, integer *, complex *, complex 
853             *, integer *, integer *), ctrsna_(char *, char *, logical *, 
854             integer *, complex *, integer *, complex *, integer *, complex *, 
855             integer *, real *, real *, integer *, integer *, complex *, 
856             integer *, real *, integer *);
857     integer minwrk, maxwrk;
858     logical wantvl, wntsnb;
859     integer hswork;
860     logical wntsne;
861     real smlnum;
862     logical lquery, wantvr, wntsnn, wntsnv;
863     extern /* Subroutine */ int ctrevc3_(char *, char *, logical *, integer *,
864              complex *, integer *, complex *, integer *, complex *, integer *,
865              integer *, integer *, complex *, integer *, real *, integer *, 
866             integer *);
867     char job[1];
868     real scl, dum[1], eps;
869     complex tmp;
870     integer lwork_trevc__;
871
872
873 /*  -- LAPACK driver routine (version 3.7.1) -- */
874 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
875 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
876 /*     June 2016 */
877
878
879 /*  ===================================================================== */
880
881
882 /*     Test the input arguments */
883
884     /* Parameter adjustments */
885     a_dim1 = *lda;
886     a_offset = 1 + a_dim1 * 1;
887     a -= a_offset;
888     --w;
889     vl_dim1 = *ldvl;
890     vl_offset = 1 + vl_dim1 * 1;
891     vl -= vl_offset;
892     vr_dim1 = *ldvr;
893     vr_offset = 1 + vr_dim1 * 1;
894     vr -= vr_offset;
895     --scale;
896     --rconde;
897     --rcondv;
898     --work;
899     --rwork;
900
901     /* Function Body */
902     *info = 0;
903     lquery = *lwork == -1;
904     wantvl = lsame_(jobvl, "V");
905     wantvr = lsame_(jobvr, "V");
906     wntsnn = lsame_(sense, "N");
907     wntsne = lsame_(sense, "E");
908     wntsnv = lsame_(sense, "V");
909     wntsnb = lsame_(sense, "B");
910     if (! (lsame_(balanc, "N") || lsame_(balanc, "S") || lsame_(balanc, "P") 
911             || lsame_(balanc, "B"))) {
912         *info = -1;
913     } else if (! wantvl && ! lsame_(jobvl, "N")) {
914         *info = -2;
915     } else if (! wantvr && ! lsame_(jobvr, "N")) {
916         *info = -3;
917     } else if (! (wntsnn || wntsne || wntsnb || wntsnv) || (wntsne || wntsnb) 
918             && ! (wantvl && wantvr)) {
919         *info = -4;
920     } else if (*n < 0) {
921         *info = -5;
922     } else if (*lda < f2cmax(1,*n)) {
923         *info = -7;
924     } else if (*ldvl < 1 || wantvl && *ldvl < *n) {
925         *info = -10;
926     } else if (*ldvr < 1 || wantvr && *ldvr < *n) {
927         *info = -12;
928     }
929
930 /*     Compute workspace */
931 /*      (Note: Comments in the code beginning "Workspace:" describe the */
932 /*       minimal amount of workspace needed at that point in the code, */
933 /*       as well as the preferred amount for good performance. */
934 /*       CWorkspace refers to complex workspace, and RWorkspace to real */
935 /*       workspace. NB refers to the optimal block size for the */
936 /*       immediately following subroutine, as returned by ILAENV. */
937 /*       HSWORK refers to the workspace preferred by CHSEQR, as */
938 /*       calculated below. HSWORK is computed assuming ILO=1 and IHI=N, */
939 /*       the worst case.) */
940
941     if (*info == 0) {
942         if (*n == 0) {
943             minwrk = 1;
944             maxwrk = 1;
945         } else {
946             maxwrk = *n + *n * ilaenv_(&c__1, "CGEHRD", " ", n, &c__1, n, &
947                     c__0, (ftnlen)6, (ftnlen)1);
948
949             if (wantvl) {
950                 ctrevc3_("L", "B", select, n, &a[a_offset], lda, &vl[
951                         vl_offset], ldvl, &vr[vr_offset], ldvr, n, &nout, &
952                         work[1], &c_n1, &rwork[1], &c_n1, &ierr);
953                 lwork_trevc__ = (integer) work[1].r;
954                 maxwrk = f2cmax(maxwrk,lwork_trevc__);
955                 chseqr_("S", "V", n, &c__1, n, &a[a_offset], lda, &w[1], &vl[
956                         vl_offset], ldvl, &work[1], &c_n1, info);
957             } else if (wantvr) {
958                 ctrevc3_("R", "B", select, n, &a[a_offset], lda, &vl[
959                         vl_offset], ldvl, &vr[vr_offset], ldvr, n, &nout, &
960                         work[1], &c_n1, &rwork[1], &c_n1, &ierr);
961                 lwork_trevc__ = (integer) work[1].r;
962                 maxwrk = f2cmax(maxwrk,lwork_trevc__);
963                 chseqr_("S", "V", n, &c__1, n, &a[a_offset], lda, &w[1], &vr[
964                         vr_offset], ldvr, &work[1], &c_n1, info);
965             } else {
966                 if (wntsnn) {
967                     chseqr_("E", "N", n, &c__1, n, &a[a_offset], lda, &w[1], &
968                             vr[vr_offset], ldvr, &work[1], &c_n1, info);
969                 } else {
970                     chseqr_("S", "N", n, &c__1, n, &a[a_offset], lda, &w[1], &
971                             vr[vr_offset], ldvr, &work[1], &c_n1, info);
972                 }
973             }
974             hswork = (integer) work[1].r;
975
976             if (! wantvl && ! wantvr) {
977                 minwrk = *n << 1;
978                 if (! (wntsnn || wntsne)) {
979 /* Computing MAX */
980                     i__1 = minwrk, i__2 = *n * *n + (*n << 1);
981                     minwrk = f2cmax(i__1,i__2);
982                 }
983                 maxwrk = f2cmax(maxwrk,hswork);
984                 if (! (wntsnn || wntsne)) {
985 /* Computing MAX */
986                     i__1 = maxwrk, i__2 = *n * *n + (*n << 1);
987                     maxwrk = f2cmax(i__1,i__2);
988                 }
989             } else {
990                 minwrk = *n << 1;
991                 if (! (wntsnn || wntsne)) {
992 /* Computing MAX */
993                     i__1 = minwrk, i__2 = *n * *n + (*n << 1);
994                     minwrk = f2cmax(i__1,i__2);
995                 }
996                 maxwrk = f2cmax(maxwrk,hswork);
997 /* Computing MAX */
998                 i__1 = maxwrk, i__2 = *n + (*n - 1) * ilaenv_(&c__1, "CUNGHR",
999                          " ", n, &c__1, n, &c_n1, (ftnlen)6, (ftnlen)1);
1000                 maxwrk = f2cmax(i__1,i__2);
1001                 if (! (wntsnn || wntsne)) {
1002 /* Computing MAX */
1003                     i__1 = maxwrk, i__2 = *n * *n + (*n << 1);
1004                     maxwrk = f2cmax(i__1,i__2);
1005                 }
1006 /* Computing MAX */
1007                 i__1 = maxwrk, i__2 = *n << 1;
1008                 maxwrk = f2cmax(i__1,i__2);
1009             }
1010             maxwrk = f2cmax(maxwrk,minwrk);
1011         }
1012         work[1].r = (real) maxwrk, work[1].i = 0.f;
1013
1014         if (*lwork < minwrk && ! lquery) {
1015             *info = -20;
1016         }
1017     }
1018
1019     if (*info != 0) {
1020         i__1 = -(*info);
1021         xerbla_("CGEEVX", &i__1, (ftnlen)6);
1022         return 0;
1023     } else if (lquery) {
1024         return 0;
1025     }
1026
1027 /*     Quick return if possible */
1028
1029     if (*n == 0) {
1030         return 0;
1031     }
1032
1033 /*     Get machine constants */
1034
1035     eps = slamch_("P");
1036     smlnum = slamch_("S");
1037     bignum = 1.f / smlnum;
1038     slabad_(&smlnum, &bignum);
1039     smlnum = sqrt(smlnum) / eps;
1040     bignum = 1.f / smlnum;
1041
1042 /*     Scale A if f2cmax element outside range [SMLNUM,BIGNUM] */
1043
1044     icond = 0;
1045     anrm = clange_("M", n, n, &a[a_offset], lda, dum);
1046     scalea = FALSE_;
1047     if (anrm > 0.f && anrm < smlnum) {
1048         scalea = TRUE_;
1049         cscale = smlnum;
1050     } else if (anrm > bignum) {
1051         scalea = TRUE_;
1052         cscale = bignum;
1053     }
1054     if (scalea) {
1055         clascl_("G", &c__0, &c__0, &anrm, &cscale, n, n, &a[a_offset], lda, &
1056                 ierr);
1057     }
1058
1059 /*     Balance the matrix and compute ABNRM */
1060
1061     cgebal_(balanc, n, &a[a_offset], lda, ilo, ihi, &scale[1], &ierr);
1062     *abnrm = clange_("1", n, n, &a[a_offset], lda, dum);
1063     if (scalea) {
1064         dum[0] = *abnrm;
1065         slascl_("G", &c__0, &c__0, &cscale, &anrm, &c__1, &c__1, dum, &c__1, &
1066                 ierr);
1067         *abnrm = dum[0];
1068     }
1069
1070 /*     Reduce to upper Hessenberg form */
1071 /*     (CWorkspace: need 2*N, prefer N+N*NB) */
1072 /*     (RWorkspace: none) */
1073
1074     itau = 1;
1075     iwrk = itau + *n;
1076     i__1 = *lwork - iwrk + 1;
1077     cgehrd_(n, ilo, ihi, &a[a_offset], lda, &work[itau], &work[iwrk], &i__1, &
1078             ierr);
1079
1080     if (wantvl) {
1081
1082 /*        Want left eigenvectors */
1083 /*        Copy Householder vectors to VL */
1084
1085         *(unsigned char *)side = 'L';
1086         clacpy_("L", n, n, &a[a_offset], lda, &vl[vl_offset], ldvl)
1087                 ;
1088
1089 /*        Generate unitary matrix in VL */
1090 /*        (CWorkspace: need 2*N-1, prefer N+(N-1)*NB) */
1091 /*        (RWorkspace: none) */
1092
1093         i__1 = *lwork - iwrk + 1;
1094         cunghr_(n, ilo, ihi, &vl[vl_offset], ldvl, &work[itau], &work[iwrk], &
1095                 i__1, &ierr);
1096
1097 /*        Perform QR iteration, accumulating Schur vectors in VL */
1098 /*        (CWorkspace: need 1, prefer HSWORK (see comments) ) */
1099 /*        (RWorkspace: none) */
1100
1101         iwrk = itau;
1102         i__1 = *lwork - iwrk + 1;
1103         chseqr_("S", "V", n, ilo, ihi, &a[a_offset], lda, &w[1], &vl[
1104                 vl_offset], ldvl, &work[iwrk], &i__1, info);
1105
1106         if (wantvr) {
1107
1108 /*           Want left and right eigenvectors */
1109 /*           Copy Schur vectors to VR */
1110
1111             *(unsigned char *)side = 'B';
1112             clacpy_("F", n, n, &vl[vl_offset], ldvl, &vr[vr_offset], ldvr);
1113         }
1114
1115     } else if (wantvr) {
1116
1117 /*        Want right eigenvectors */
1118 /*        Copy Householder vectors to VR */
1119
1120         *(unsigned char *)side = 'R';
1121         clacpy_("L", n, n, &a[a_offset], lda, &vr[vr_offset], ldvr)
1122                 ;
1123
1124 /*        Generate unitary matrix in VR */
1125 /*        (CWorkspace: need 2*N-1, prefer N+(N-1)*NB) */
1126 /*        (RWorkspace: none) */
1127
1128         i__1 = *lwork - iwrk + 1;
1129         cunghr_(n, ilo, ihi, &vr[vr_offset], ldvr, &work[itau], &work[iwrk], &
1130                 i__1, &ierr);
1131
1132 /*        Perform QR iteration, accumulating Schur vectors in VR */
1133 /*        (CWorkspace: need 1, prefer HSWORK (see comments) ) */
1134 /*        (RWorkspace: none) */
1135
1136         iwrk = itau;
1137         i__1 = *lwork - iwrk + 1;
1138         chseqr_("S", "V", n, ilo, ihi, &a[a_offset], lda, &w[1], &vr[
1139                 vr_offset], ldvr, &work[iwrk], &i__1, info);
1140
1141     } else {
1142
1143 /*        Compute eigenvalues only */
1144 /*        If condition numbers desired, compute Schur form */
1145
1146         if (wntsnn) {
1147             *(unsigned char *)job = 'E';
1148         } else {
1149             *(unsigned char *)job = 'S';
1150         }
1151
1152 /*        (CWorkspace: need 1, prefer HSWORK (see comments) ) */
1153 /*        (RWorkspace: none) */
1154
1155         iwrk = itau;
1156         i__1 = *lwork - iwrk + 1;
1157         chseqr_(job, "N", n, ilo, ihi, &a[a_offset], lda, &w[1], &vr[
1158                 vr_offset], ldvr, &work[iwrk], &i__1, info);
1159     }
1160
1161 /*     If INFO .NE. 0 from CHSEQR, then quit */
1162
1163     if (*info != 0) {
1164         goto L50;
1165     }
1166
1167     if (wantvl || wantvr) {
1168
1169 /*        Compute left and/or right eigenvectors */
1170 /*        (CWorkspace: need 2*N, prefer N + 2*N*NB) */
1171 /*        (RWorkspace: need N) */
1172
1173         i__1 = *lwork - iwrk + 1;
1174         ctrevc3_(side, "B", select, n, &a[a_offset], lda, &vl[vl_offset], 
1175                 ldvl, &vr[vr_offset], ldvr, n, &nout, &work[iwrk], &i__1, &
1176                 rwork[1], n, &ierr);
1177     }
1178
1179 /*     Compute condition numbers if desired */
1180 /*     (CWorkspace: need N*N+2*N unless SENSE = 'E') */
1181 /*     (RWorkspace: need 2*N unless SENSE = 'E') */
1182
1183     if (! wntsnn) {
1184         ctrsna_(sense, "A", select, n, &a[a_offset], lda, &vl[vl_offset], 
1185                 ldvl, &vr[vr_offset], ldvr, &rconde[1], &rcondv[1], n, &nout, 
1186                 &work[iwrk], n, &rwork[1], &icond);
1187     }
1188
1189     if (wantvl) {
1190
1191 /*        Undo balancing of left eigenvectors */
1192
1193         cgebak_(balanc, "L", n, ilo, ihi, &scale[1], n, &vl[vl_offset], ldvl, 
1194                 &ierr);
1195
1196 /*        Normalize left eigenvectors and make largest component real */
1197
1198         i__1 = *n;
1199         for (i__ = 1; i__ <= i__1; ++i__) {
1200             scl = 1.f / scnrm2_(n, &vl[i__ * vl_dim1 + 1], &c__1);
1201             csscal_(n, &scl, &vl[i__ * vl_dim1 + 1], &c__1);
1202             i__2 = *n;
1203             for (k = 1; k <= i__2; ++k) {
1204                 i__3 = k + i__ * vl_dim1;
1205 /* Computing 2nd power */
1206                 r__1 = vl[i__3].r;
1207 /* Computing 2nd power */
1208                 r__2 = r_imag(&vl[k + i__ * vl_dim1]);
1209                 rwork[k] = r__1 * r__1 + r__2 * r__2;
1210 /* L10: */
1211             }
1212             k = isamax_(n, &rwork[1], &c__1);
1213             r_cnjg(&q__2, &vl[k + i__ * vl_dim1]);
1214             r__1 = sqrt(rwork[k]);
1215             q__1.r = q__2.r / r__1, q__1.i = q__2.i / r__1;
1216             tmp.r = q__1.r, tmp.i = q__1.i;
1217             cscal_(n, &tmp, &vl[i__ * vl_dim1 + 1], &c__1);
1218             i__2 = k + i__ * vl_dim1;
1219             i__3 = k + i__ * vl_dim1;
1220             r__1 = vl[i__3].r;
1221             q__1.r = r__1, q__1.i = 0.f;
1222             vl[i__2].r = q__1.r, vl[i__2].i = q__1.i;
1223 /* L20: */
1224         }
1225     }
1226
1227     if (wantvr) {
1228
1229 /*        Undo balancing of right eigenvectors */
1230
1231         cgebak_(balanc, "R", n, ilo, ihi, &scale[1], n, &vr[vr_offset], ldvr, 
1232                 &ierr);
1233
1234 /*        Normalize right eigenvectors and make largest component real */
1235
1236         i__1 = *n;
1237         for (i__ = 1; i__ <= i__1; ++i__) {
1238             scl = 1.f / scnrm2_(n, &vr[i__ * vr_dim1 + 1], &c__1);
1239             csscal_(n, &scl, &vr[i__ * vr_dim1 + 1], &c__1);
1240             i__2 = *n;
1241             for (k = 1; k <= i__2; ++k) {
1242                 i__3 = k + i__ * vr_dim1;
1243 /* Computing 2nd power */
1244                 r__1 = vr[i__3].r;
1245 /* Computing 2nd power */
1246                 r__2 = r_imag(&vr[k + i__ * vr_dim1]);
1247                 rwork[k] = r__1 * r__1 + r__2 * r__2;
1248 /* L30: */
1249             }
1250             k = isamax_(n, &rwork[1], &c__1);
1251             r_cnjg(&q__2, &vr[k + i__ * vr_dim1]);
1252             r__1 = sqrt(rwork[k]);
1253             q__1.r = q__2.r / r__1, q__1.i = q__2.i / r__1;
1254             tmp.r = q__1.r, tmp.i = q__1.i;
1255             cscal_(n, &tmp, &vr[i__ * vr_dim1 + 1], &c__1);
1256             i__2 = k + i__ * vr_dim1;
1257             i__3 = k + i__ * vr_dim1;
1258             r__1 = vr[i__3].r;
1259             q__1.r = r__1, q__1.i = 0.f;
1260             vr[i__2].r = q__1.r, vr[i__2].i = q__1.i;
1261 /* L40: */
1262         }
1263     }
1264
1265 /*     Undo scaling if necessary */
1266
1267 L50:
1268     if (scalea) {
1269         i__1 = *n - *info;
1270 /* Computing MAX */
1271         i__3 = *n - *info;
1272         i__2 = f2cmax(i__3,1);
1273         clascl_("G", &c__0, &c__0, &cscale, &anrm, &i__1, &c__1, &w[*info + 1]
1274                 , &i__2, &ierr);
1275         if (*info == 0) {
1276             if ((wntsnv || wntsnb) && icond == 0) {
1277                 slascl_("G", &c__0, &c__0, &cscale, &anrm, n, &c__1, &rcondv[
1278                         1], n, &ierr);
1279             }
1280         } else {
1281             i__1 = *ilo - 1;
1282             clascl_("G", &c__0, &c__0, &cscale, &anrm, &i__1, &c__1, &w[1], n,
1283                      &ierr);
1284         }
1285     }
1286
1287     work[1].r = (real) maxwrk, work[1].i = 0.f;
1288     return 0;
1289
1290 /*     End of CGEEVX */
1291
1292 } /* cgeevx_ */
1293