C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / DEPRECATED / sgegs.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_n1 = -1;
517 static real c_b36 = 0.f;
518 static real c_b37 = 1.f;
519
520 /* > \brief <b> SGEGS computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matr
521 ices</b> */
522
523 /*  =========== DOCUMENTATION =========== */
524
525 /* Online html documentation available at */
526 /*            http://www.netlib.org/lapack/explore-html/ */
527
528 /* > \htmlonly */
529 /* > Download SGEGS + dependencies */
530 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgegs.f
531 "> */
532 /* > [TGZ]</a> */
533 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgegs.f
534 "> */
535 /* > [ZIP]</a> */
536 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgegs.f
537 "> */
538 /* > [TXT]</a> */
539 /* > \endhtmlonly */
540
541 /*  Definition: */
542 /*  =========== */
543
544 /*       SUBROUTINE SGEGS( JOBVSL, JOBVSR, N, A, LDA, B, LDB, ALPHAR, */
545 /*                         ALPHAI, BETA, VSL, LDVSL, VSR, LDVSR, WORK, */
546 /*                         LWORK, INFO ) */
547
548 /*       CHARACTER          JOBVSL, JOBVSR */
549 /*       INTEGER            INFO, LDA, LDB, LDVSL, LDVSR, LWORK, N */
550 /*       REAL               A( LDA, * ), ALPHAI( * ), ALPHAR( * ), */
551 /*      $                   B( LDB, * ), BETA( * ), VSL( LDVSL, * ), */
552 /*      $                   VSR( LDVSR, * ), WORK( * ) */
553
554
555 /* > \par Purpose: */
556 /*  ============= */
557 /* > */
558 /* > \verbatim */
559 /* > */
560 /* > This routine is deprecated and has been replaced by routine SGGES. */
561 /* > */
562 /* > SGEGS computes the eigenvalues, real Schur form, and, optionally, */
563 /* > left and or/right Schur vectors of a real matrix pair (A,B). */
564 /* > Given two square matrices A and B, the generalized real Schur */
565 /* > factorization has the form */
566 /* > */
567 /* >   A = Q*S*Z**T,  B = Q*T*Z**T */
568 /* > */
569 /* > where Q and Z are orthogonal matrices, T is upper triangular, and S */
570 /* > is an upper quasi-triangular matrix with 1-by-1 and 2-by-2 diagonal */
571 /* > blocks, the 2-by-2 blocks corresponding to complex conjugate pairs */
572 /* > of eigenvalues of (A,B).  The columns of Q are the left Schur vectors */
573 /* > and the columns of Z are the right Schur vectors. */
574 /* > */
575 /* > If only the eigenvalues of (A,B) are needed, the driver routine */
576 /* > SGEGV should be used instead.  See SGEGV for a description of the */
577 /* > eigenvalues of the generalized nonsymmetric eigenvalue problem */
578 /* > (GNEP). */
579 /* > \endverbatim */
580
581 /*  Arguments: */
582 /*  ========== */
583
584 /* > \param[in] JOBVSL */
585 /* > \verbatim */
586 /* >          JOBVSL is CHARACTER*1 */
587 /* >          = 'N':  do not compute the left Schur vectors; */
588 /* >          = 'V':  compute the left Schur vectors (returned in VSL). */
589 /* > \endverbatim */
590 /* > */
591 /* > \param[in] JOBVSR */
592 /* > \verbatim */
593 /* >          JOBVSR is CHARACTER*1 */
594 /* >          = 'N':  do not compute the right Schur vectors; */
595 /* >          = 'V':  compute the right Schur vectors (returned in VSR). */
596 /* > \endverbatim */
597 /* > */
598 /* > \param[in] N */
599 /* > \verbatim */
600 /* >          N is INTEGER */
601 /* >          The order of the matrices A, B, VSL, and VSR.  N >= 0. */
602 /* > \endverbatim */
603 /* > */
604 /* > \param[in,out] A */
605 /* > \verbatim */
606 /* >          A is REAL array, dimension (LDA, N) */
607 /* >          On entry, the matrix A. */
608 /* >          On exit, the upper quasi-triangular matrix S from the */
609 /* >          generalized real Schur factorization. */
610 /* > \endverbatim */
611 /* > */
612 /* > \param[in] LDA */
613 /* > \verbatim */
614 /* >          LDA is INTEGER */
615 /* >          The leading dimension of A.  LDA >= f2cmax(1,N). */
616 /* > \endverbatim */
617 /* > */
618 /* > \param[in,out] B */
619 /* > \verbatim */
620 /* >          B is REAL array, dimension (LDB, N) */
621 /* >          On entry, the matrix B. */
622 /* >          On exit, the upper triangular matrix T from the generalized */
623 /* >          real Schur factorization. */
624 /* > \endverbatim */
625 /* > */
626 /* > \param[in] LDB */
627 /* > \verbatim */
628 /* >          LDB is INTEGER */
629 /* >          The leading dimension of B.  LDB >= f2cmax(1,N). */
630 /* > \endverbatim */
631 /* > */
632 /* > \param[out] ALPHAR */
633 /* > \verbatim */
634 /* >          ALPHAR is REAL array, dimension (N) */
635 /* >          The real parts of each scalar alpha defining an eigenvalue */
636 /* >          of GNEP. */
637 /* > \endverbatim */
638 /* > */
639 /* > \param[out] ALPHAI */
640 /* > \verbatim */
641 /* >          ALPHAI is REAL array, dimension (N) */
642 /* >          The imaginary parts of each scalar alpha defining an */
643 /* >          eigenvalue of GNEP.  If ALPHAI(j) is zero, then the j-th */
644 /* >          eigenvalue is real; if positive, then the j-th and (j+1)-st */
645 /* >          eigenvalues are a complex conjugate pair, with */
646 /* >          ALPHAI(j+1) = -ALPHAI(j). */
647 /* > \endverbatim */
648 /* > */
649 /* > \param[out] BETA */
650 /* > \verbatim */
651 /* >          BETA is REAL array, dimension (N) */
652 /* >          The scalars beta that define the eigenvalues of GNEP. */
653 /* >          Together, the quantities alpha = (ALPHAR(j),ALPHAI(j)) and */
654 /* >          beta = BETA(j) represent the j-th eigenvalue of the matrix */
655 /* >          pair (A,B), in one of the forms lambda = alpha/beta or */
656 /* >          mu = beta/alpha.  Since either lambda or mu may overflow, */
657 /* >          they should not, in general, be computed. */
658 /* > \endverbatim */
659 /* > */
660 /* > \param[out] VSL */
661 /* > \verbatim */
662 /* >          VSL is REAL array, dimension (LDVSL,N) */
663 /* >          If JOBVSL = 'V', the matrix of left Schur vectors Q. */
664 /* >          Not referenced if JOBVSL = 'N'. */
665 /* > \endverbatim */
666 /* > */
667 /* > \param[in] LDVSL */
668 /* > \verbatim */
669 /* >          LDVSL is INTEGER */
670 /* >          The leading dimension of the matrix VSL. LDVSL >=1, and */
671 /* >          if JOBVSL = 'V', LDVSL >= N. */
672 /* > \endverbatim */
673 /* > */
674 /* > \param[out] VSR */
675 /* > \verbatim */
676 /* >          VSR is REAL array, dimension (LDVSR,N) */
677 /* >          If JOBVSR = 'V', the matrix of right Schur vectors Z. */
678 /* >          Not referenced if JOBVSR = 'N'. */
679 /* > \endverbatim */
680 /* > */
681 /* > \param[in] LDVSR */
682 /* > \verbatim */
683 /* >          LDVSR is INTEGER */
684 /* >          The leading dimension of the matrix VSR. LDVSR >= 1, and */
685 /* >          if JOBVSR = 'V', LDVSR >= N. */
686 /* > \endverbatim */
687 /* > */
688 /* > \param[out] WORK */
689 /* > \verbatim */
690 /* >          WORK is REAL array, dimension (MAX(1,LWORK)) */
691 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
692 /* > \endverbatim */
693 /* > */
694 /* > \param[in] LWORK */
695 /* > \verbatim */
696 /* >          LWORK is INTEGER */
697 /* >          The dimension of the array WORK.  LWORK >= f2cmax(1,4*N). */
698 /* >          For good performance, LWORK must generally be larger. */
699 /* >          To compute the optimal value of LWORK, call ILAENV to get */
700 /* >          blocksizes (for SGEQRF, SORMQR, and SORGQR.)  Then compute: */
701 /* >          NB  -- MAX of the blocksizes for SGEQRF, SORMQR, and SORGQR */
702 /* >          The optimal LWORK is  2*N + N*(NB+1). */
703 /* > */
704 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
705 /* >          only calculates the optimal size of the WORK array, returns */
706 /* >          this value as the first entry of the WORK array, and no error */
707 /* >          message related to LWORK is issued by XERBLA. */
708 /* > \endverbatim */
709 /* > */
710 /* > \param[out] INFO */
711 /* > \verbatim */
712 /* >          INFO is INTEGER */
713 /* >          = 0:  successful exit */
714 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value. */
715 /* >          = 1,...,N: */
716 /* >                The QZ iteration failed.  (A,B) are not in Schur */
717 /* >                form, but ALPHAR(j), ALPHAI(j), and BETA(j) should */
718 /* >                be correct for j=INFO+1,...,N. */
719 /* >          > N:  errors that usually indicate LAPACK problems: */
720 /* >                =N+1: error return from SGGBAL */
721 /* >                =N+2: error return from SGEQRF */
722 /* >                =N+3: error return from SORMQR */
723 /* >                =N+4: error return from SORGQR */
724 /* >                =N+5: error return from SGGHRD */
725 /* >                =N+6: error return from SHGEQZ (other than failed */
726 /* >                                                iteration) */
727 /* >                =N+7: error return from SGGBAK (computing VSL) */
728 /* >                =N+8: error return from SGGBAK (computing VSR) */
729 /* >                =N+9: error return from SLASCL (various places) */
730 /* > \endverbatim */
731
732 /*  Authors: */
733 /*  ======== */
734
735 /* > \author Univ. of Tennessee */
736 /* > \author Univ. of California Berkeley */
737 /* > \author Univ. of Colorado Denver */
738 /* > \author NAG Ltd. */
739
740 /* > \date December 2016 */
741
742 /* > \ingroup realGEeigen */
743
744 /*  ===================================================================== */
745 /* Subroutine */ int sgegs_(char *jobvsl, char *jobvsr, integer *n, real *a, 
746         integer *lda, real *b, integer *ldb, real *alphar, real *alphai, real 
747         *beta, real *vsl, integer *ldvsl, real *vsr, integer *ldvsr, real *
748         work, integer *lwork, integer *info)
749 {
750     /* System generated locals */
751     integer a_dim1, a_offset, b_dim1, b_offset, vsl_dim1, vsl_offset, 
752             vsr_dim1, vsr_offset, i__1, i__2;
753
754     /* Local variables */
755     real anrm, bnrm;
756     integer itau, lopt;
757     extern logical lsame_(char *, char *);
758     integer ileft, iinfo, icols;
759     logical ilvsl;
760     integer iwork;
761     logical ilvsr;
762     integer irows, nb;
763     extern /* Subroutine */ int sggbak_(char *, char *, integer *, integer *, 
764             integer *, real *, real *, integer *, real *, integer *, integer *
765             ), sggbal_(char *, integer *, real *, integer *, 
766             real *, integer *, integer *, integer *, real *, real *, real *, 
767             integer *);
768     logical ilascl, ilbscl;
769     extern real slamch_(char *), slange_(char *, integer *, integer *,
770              real *, integer *, real *);
771     real safmin;
772     extern /* Subroutine */ int sgghrd_(char *, char *, integer *, integer *, 
773             integer *, real *, integer *, real *, integer *, real *, integer *
774             , real *, integer *, integer *), xerbla_(char *, 
775             integer *);
776     extern integer ilaenv_(integer *, char *, char *, integer *, integer *, 
777             integer *, integer *, ftnlen, ftnlen);
778     real bignum;
779     extern /* Subroutine */ int slascl_(char *, integer *, integer *, real *, 
780             real *, integer *, integer *, real *, integer *, integer *);
781     integer ijobvl, iright;
782     extern /* Subroutine */ int sgeqrf_(integer *, integer *, real *, integer 
783             *, real *, real *, integer *, integer *);
784     integer ijobvr;
785     extern /* Subroutine */ int slacpy_(char *, integer *, integer *, real *, 
786             integer *, real *, integer *), slaset_(char *, integer *, 
787             integer *, real *, real *, real *, integer *);
788     real anrmto;
789     integer lwkmin, nb1, nb2, nb3;
790     real bnrmto;
791     extern /* Subroutine */ int shgeqz_(char *, char *, char *, integer *, 
792             integer *, integer *, real *, integer *, real *, integer *, real *
793             , real *, real *, real *, integer *, real *, integer *, real *, 
794             integer *, integer *);
795     real smlnum;
796     extern /* Subroutine */ int sorgqr_(integer *, integer *, integer *, real 
797             *, integer *, real *, real *, integer *, integer *);
798     integer lwkopt;
799     logical lquery;
800     extern /* Subroutine */ int sormqr_(char *, char *, integer *, integer *, 
801             integer *, real *, integer *, real *, real *, integer *, real *, 
802             integer *, integer *);
803     integer ihi, ilo;
804     real eps;
805
806
807 /*  -- LAPACK driver routine (version 3.7.0) -- */
808 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
809 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
810 /*     December 2016 */
811
812
813 /*  ===================================================================== */
814
815
816 /*     Decode the input arguments */
817
818     /* Parameter adjustments */
819     a_dim1 = *lda;
820     a_offset = 1 + a_dim1 * 1;
821     a -= a_offset;
822     b_dim1 = *ldb;
823     b_offset = 1 + b_dim1 * 1;
824     b -= b_offset;
825     --alphar;
826     --alphai;
827     --beta;
828     vsl_dim1 = *ldvsl;
829     vsl_offset = 1 + vsl_dim1 * 1;
830     vsl -= vsl_offset;
831     vsr_dim1 = *ldvsr;
832     vsr_offset = 1 + vsr_dim1 * 1;
833     vsr -= vsr_offset;
834     --work;
835
836     /* Function Body */
837     if (lsame_(jobvsl, "N")) {
838         ijobvl = 1;
839         ilvsl = FALSE_;
840     } else if (lsame_(jobvsl, "V")) {
841         ijobvl = 2;
842         ilvsl = TRUE_;
843     } else {
844         ijobvl = -1;
845         ilvsl = FALSE_;
846     }
847
848     if (lsame_(jobvsr, "N")) {
849         ijobvr = 1;
850         ilvsr = FALSE_;
851     } else if (lsame_(jobvsr, "V")) {
852         ijobvr = 2;
853         ilvsr = TRUE_;
854     } else {
855         ijobvr = -1;
856         ilvsr = FALSE_;
857     }
858
859 /*     Test the input arguments */
860
861 /* Computing MAX */
862     i__1 = *n << 2;
863     lwkmin = f2cmax(i__1,1);
864     lwkopt = lwkmin;
865     work[1] = (real) lwkopt;
866     lquery = *lwork == -1;
867     *info = 0;
868     if (ijobvl <= 0) {
869         *info = -1;
870     } else if (ijobvr <= 0) {
871         *info = -2;
872     } else if (*n < 0) {
873         *info = -3;
874     } else if (*lda < f2cmax(1,*n)) {
875         *info = -5;
876     } else if (*ldb < f2cmax(1,*n)) {
877         *info = -7;
878     } else if (*ldvsl < 1 || ilvsl && *ldvsl < *n) {
879         *info = -12;
880     } else if (*ldvsr < 1 || ilvsr && *ldvsr < *n) {
881         *info = -14;
882     } else if (*lwork < lwkmin && ! lquery) {
883         *info = -16;
884     }
885
886     if (*info == 0) {
887         nb1 = ilaenv_(&c__1, "SGEQRF", " ", n, n, &c_n1, &c_n1, (ftnlen)6, (
888                 ftnlen)1);
889         nb2 = ilaenv_(&c__1, "SORMQR", " ", n, n, n, &c_n1, (ftnlen)6, (
890                 ftnlen)1);
891         nb3 = ilaenv_(&c__1, "SORGQR", " ", n, n, n, &c_n1, (ftnlen)6, (
892                 ftnlen)1);
893 /* Computing MAX */
894         i__1 = f2cmax(nb1,nb2);
895         nb = f2cmax(i__1,nb3);
896         lopt = (*n << 1) + *n * (nb + 1);
897         work[1] = (real) lopt;
898     }
899
900     if (*info != 0) {
901         i__1 = -(*info);
902         xerbla_("SGEGS ", &i__1);
903         return 0;
904     } else if (lquery) {
905         return 0;
906     }
907
908 /*     Quick return if possible */
909
910     if (*n == 0) {
911         return 0;
912     }
913
914 /*     Get machine constants */
915
916     eps = slamch_("E") * slamch_("B");
917     safmin = slamch_("S");
918     smlnum = *n * safmin / eps;
919     bignum = 1.f / smlnum;
920
921 /*     Scale A if f2cmax element outside range [SMLNUM,BIGNUM] */
922
923     anrm = slange_("M", n, n, &a[a_offset], lda, &work[1]);
924     ilascl = FALSE_;
925     if (anrm > 0.f && anrm < smlnum) {
926         anrmto = smlnum;
927         ilascl = TRUE_;
928     } else if (anrm > bignum) {
929         anrmto = bignum;
930         ilascl = TRUE_;
931     }
932
933     if (ilascl) {
934         slascl_("G", &c_n1, &c_n1, &anrm, &anrmto, n, n, &a[a_offset], lda, &
935                 iinfo);
936         if (iinfo != 0) {
937             *info = *n + 9;
938             return 0;
939         }
940     }
941
942 /*     Scale B if f2cmax element outside range [SMLNUM,BIGNUM] */
943
944     bnrm = slange_("M", n, n, &b[b_offset], ldb, &work[1]);
945     ilbscl = FALSE_;
946     if (bnrm > 0.f && bnrm < smlnum) {
947         bnrmto = smlnum;
948         ilbscl = TRUE_;
949     } else if (bnrm > bignum) {
950         bnrmto = bignum;
951         ilbscl = TRUE_;
952     }
953
954     if (ilbscl) {
955         slascl_("G", &c_n1, &c_n1, &bnrm, &bnrmto, n, n, &b[b_offset], ldb, &
956                 iinfo);
957         if (iinfo != 0) {
958             *info = *n + 9;
959             return 0;
960         }
961     }
962
963 /*     Permute the matrix to make it more nearly triangular */
964 /*     Workspace layout:  (2*N words -- "work..." not actually used) */
965 /*        left_permutation, right_permutation, work... */
966
967     ileft = 1;
968     iright = *n + 1;
969     iwork = iright + *n;
970     sggbal_("P", n, &a[a_offset], lda, &b[b_offset], ldb, &ilo, &ihi, &work[
971             ileft], &work[iright], &work[iwork], &iinfo);
972     if (iinfo != 0) {
973         *info = *n + 1;
974         goto L10;
975     }
976
977 /*     Reduce B to triangular form, and initialize VSL and/or VSR */
978 /*     Workspace layout:  ("work..." must have at least N words) */
979 /*        left_permutation, right_permutation, tau, work... */
980
981     irows = ihi + 1 - ilo;
982     icols = *n + 1 - ilo;
983     itau = iwork;
984     iwork = itau + irows;
985     i__1 = *lwork + 1 - iwork;
986     sgeqrf_(&irows, &icols, &b[ilo + ilo * b_dim1], ldb, &work[itau], &work[
987             iwork], &i__1, &iinfo);
988     if (iinfo >= 0) {
989 /* Computing MAX */
990         i__1 = lwkopt, i__2 = (integer) work[iwork] + iwork - 1;
991         lwkopt = f2cmax(i__1,i__2);
992     }
993     if (iinfo != 0) {
994         *info = *n + 2;
995         goto L10;
996     }
997
998     i__1 = *lwork + 1 - iwork;
999     sormqr_("L", "T", &irows, &icols, &irows, &b[ilo + ilo * b_dim1], ldb, &
1000             work[itau], &a[ilo + ilo * a_dim1], lda, &work[iwork], &i__1, &
1001             iinfo);
1002     if (iinfo >= 0) {
1003 /* Computing MAX */
1004         i__1 = lwkopt, i__2 = (integer) work[iwork] + iwork - 1;
1005         lwkopt = f2cmax(i__1,i__2);
1006     }
1007     if (iinfo != 0) {
1008         *info = *n + 3;
1009         goto L10;
1010     }
1011
1012     if (ilvsl) {
1013         slaset_("Full", n, n, &c_b36, &c_b37, &vsl[vsl_offset], ldvsl);
1014         i__1 = irows - 1;
1015         i__2 = irows - 1;
1016         slacpy_("L", &i__1, &i__2, &b[ilo + 1 + ilo * b_dim1], ldb, &vsl[ilo 
1017                 + 1 + ilo * vsl_dim1], ldvsl);
1018         i__1 = *lwork + 1 - iwork;
1019         sorgqr_(&irows, &irows, &irows, &vsl[ilo + ilo * vsl_dim1], ldvsl, &
1020                 work[itau], &work[iwork], &i__1, &iinfo);
1021         if (iinfo >= 0) {
1022 /* Computing MAX */
1023             i__1 = lwkopt, i__2 = (integer) work[iwork] + iwork - 1;
1024             lwkopt = f2cmax(i__1,i__2);
1025         }
1026         if (iinfo != 0) {
1027             *info = *n + 4;
1028             goto L10;
1029         }
1030     }
1031
1032     if (ilvsr) {
1033         slaset_("Full", n, n, &c_b36, &c_b37, &vsr[vsr_offset], ldvsr);
1034     }
1035
1036 /*     Reduce to generalized Hessenberg form */
1037
1038     sgghrd_(jobvsl, jobvsr, n, &ilo, &ihi, &a[a_offset], lda, &b[b_offset], 
1039             ldb, &vsl[vsl_offset], ldvsl, &vsr[vsr_offset], ldvsr, &iinfo);
1040     if (iinfo != 0) {
1041         *info = *n + 5;
1042         goto L10;
1043     }
1044
1045 /*     Perform QZ algorithm, computing Schur vectors if desired */
1046 /*     Workspace layout:  ("work..." must have at least 1 word) */
1047 /*        left_permutation, right_permutation, work... */
1048
1049     iwork = itau;
1050     i__1 = *lwork + 1 - iwork;
1051     shgeqz_("S", jobvsl, jobvsr, n, &ilo, &ihi, &a[a_offset], lda, &b[
1052             b_offset], ldb, &alphar[1], &alphai[1], &beta[1], &vsl[vsl_offset]
1053             , ldvsl, &vsr[vsr_offset], ldvsr, &work[iwork], &i__1, &iinfo);
1054     if (iinfo >= 0) {
1055 /* Computing MAX */
1056         i__1 = lwkopt, i__2 = (integer) work[iwork] + iwork - 1;
1057         lwkopt = f2cmax(i__1,i__2);
1058     }
1059     if (iinfo != 0) {
1060         if (iinfo > 0 && iinfo <= *n) {
1061             *info = iinfo;
1062         } else if (iinfo > *n && iinfo <= *n << 1) {
1063             *info = iinfo - *n;
1064         } else {
1065             *info = *n + 6;
1066         }
1067         goto L10;
1068     }
1069
1070 /*     Apply permutation to VSL and VSR */
1071
1072     if (ilvsl) {
1073         sggbak_("P", "L", n, &ilo, &ihi, &work[ileft], &work[iright], n, &vsl[
1074                 vsl_offset], ldvsl, &iinfo);
1075         if (iinfo != 0) {
1076             *info = *n + 7;
1077             goto L10;
1078         }
1079     }
1080     if (ilvsr) {
1081         sggbak_("P", "R", n, &ilo, &ihi, &work[ileft], &work[iright], n, &vsr[
1082                 vsr_offset], ldvsr, &iinfo);
1083         if (iinfo != 0) {
1084             *info = *n + 8;
1085             goto L10;
1086         }
1087     }
1088
1089 /*     Undo scaling */
1090
1091     if (ilascl) {
1092         slascl_("H", &c_n1, &c_n1, &anrmto, &anrm, n, n, &a[a_offset], lda, &
1093                 iinfo);
1094         if (iinfo != 0) {
1095             *info = *n + 9;
1096             return 0;
1097         }
1098         slascl_("G", &c_n1, &c_n1, &anrmto, &anrm, n, &c__1, &alphar[1], n, &
1099                 iinfo);
1100         if (iinfo != 0) {
1101             *info = *n + 9;
1102             return 0;
1103         }
1104         slascl_("G", &c_n1, &c_n1, &anrmto, &anrm, n, &c__1, &alphai[1], n, &
1105                 iinfo);
1106         if (iinfo != 0) {
1107             *info = *n + 9;
1108             return 0;
1109         }
1110     }
1111
1112     if (ilbscl) {
1113         slascl_("U", &c_n1, &c_n1, &bnrmto, &bnrm, n, n, &b[b_offset], ldb, &
1114                 iinfo);
1115         if (iinfo != 0) {
1116             *info = *n + 9;
1117             return 0;
1118         }
1119         slascl_("G", &c_n1, &c_n1, &bnrmto, &bnrm, n, &c__1, &beta[1], n, &
1120                 iinfo);
1121         if (iinfo != 0) {
1122             *info = *n + 9;
1123             return 0;
1124         }
1125     }
1126
1127 L10:
1128     work[1] = (real) lwkopt;
1129
1130     return 0;
1131
1132 /*     End of SGEGS */
1133
1134 } /* sgegs_ */
1135