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