C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / chegvx.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 complex c_b1 = {1.f,0.f};
516 static integer c__1 = 1;
517 static integer c_n1 = -1;
518
519 /* > \brief \b CHEGVX */
520
521 /*  =========== DOCUMENTATION =========== */
522
523 /* Online html documentation available at */
524 /*            http://www.netlib.org/lapack/explore-html/ */
525
526 /* > \htmlonly */
527 /* > Download CHEGVX + dependencies */
528 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/chegvx.
529 f"> */
530 /* > [TGZ]</a> */
531 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/chegvx.
532 f"> */
533 /* > [ZIP]</a> */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/chegvx.
535 f"> */
536 /* > [TXT]</a> */
537 /* > \endhtmlonly */
538
539 /*  Definition: */
540 /*  =========== */
541
542 /*       SUBROUTINE CHEGVX( ITYPE, JOBZ, RANGE, UPLO, N, A, LDA, B, LDB, */
543 /*                          VL, VU, IL, IU, ABSTOL, M, W, Z, LDZ, WORK, */
544 /*                          LWORK, RWORK, IWORK, IFAIL, INFO ) */
545
546 /*       CHARACTER          JOBZ, RANGE, UPLO */
547 /*       INTEGER            IL, INFO, ITYPE, IU, LDA, LDB, LDZ, LWORK, M, N */
548 /*       REAL               ABSTOL, VL, VU */
549 /*       INTEGER            IFAIL( * ), IWORK( * ) */
550 /*       REAL               RWORK( * ), W( * ) */
551 /*       COMPLEX            A( LDA, * ), B( LDB, * ), WORK( * ), */
552 /*      $                   Z( LDZ, * ) */
553
554
555 /* > \par Purpose: */
556 /*  ============= */
557 /* > */
558 /* > \verbatim */
559 /* > */
560 /* > CHEGVX computes selected eigenvalues, and optionally, eigenvectors */
561 /* > of a complex generalized Hermitian-definite eigenproblem, of the form */
562 /* > A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.  Here A and */
563 /* > B are assumed to be Hermitian and B is also positive definite. */
564 /* > Eigenvalues and eigenvectors can be selected by specifying either a */
565 /* > range of values or a range of indices for the desired eigenvalues. */
566 /* > \endverbatim */
567
568 /*  Arguments: */
569 /*  ========== */
570
571 /* > \param[in] ITYPE */
572 /* > \verbatim */
573 /* >          ITYPE is INTEGER */
574 /* >          Specifies the problem type to be solved: */
575 /* >          = 1:  A*x = (lambda)*B*x */
576 /* >          = 2:  A*B*x = (lambda)*x */
577 /* >          = 3:  B*A*x = (lambda)*x */
578 /* > \endverbatim */
579 /* > */
580 /* > \param[in] JOBZ */
581 /* > \verbatim */
582 /* >          JOBZ is CHARACTER*1 */
583 /* >          = 'N':  Compute eigenvalues only; */
584 /* >          = 'V':  Compute eigenvalues and eigenvectors. */
585 /* > \endverbatim */
586 /* > */
587 /* > \param[in] RANGE */
588 /* > \verbatim */
589 /* >          RANGE is CHARACTER*1 */
590 /* >          = 'A': all eigenvalues will be found. */
591 /* >          = 'V': all eigenvalues in the half-open interval (VL,VU] */
592 /* >                 will be found. */
593 /* >          = 'I': the IL-th through IU-th eigenvalues will be found. */
594 /* > \endverbatim */
595 /* > */
596 /* > \param[in] UPLO */
597 /* > \verbatim */
598 /* >          UPLO is CHARACTER*1 */
599 /* >          = 'U':  Upper triangles of A and B are stored; */
600 /* >          = 'L':  Lower triangles of A and B are stored. */
601 /* > \endverbatim */
602 /* > */
603 /* > \param[in] N */
604 /* > \verbatim */
605 /* >          N is INTEGER */
606 /* >          The order of the matrices A and B.  N >= 0. */
607 /* > \endverbatim */
608 /* > */
609 /* > \param[in,out] A */
610 /* > \verbatim */
611 /* >          A is COMPLEX array, dimension (LDA, N) */
612 /* >          On entry, the Hermitian matrix A.  If UPLO = 'U', the */
613 /* >          leading N-by-N upper triangular part of A contains the */
614 /* >          upper triangular part of the matrix A.  If UPLO = 'L', */
615 /* >          the leading N-by-N lower triangular part of A contains */
616 /* >          the lower triangular part of the matrix A. */
617 /* > */
618 /* >          On exit,  the lower triangle (if UPLO='L') or the upper */
619 /* >          triangle (if UPLO='U') of A, including the diagonal, is */
620 /* >          destroyed. */
621 /* > \endverbatim */
622 /* > */
623 /* > \param[in] LDA */
624 /* > \verbatim */
625 /* >          LDA is INTEGER */
626 /* >          The leading dimension of the array A.  LDA >= f2cmax(1,N). */
627 /* > \endverbatim */
628 /* > */
629 /* > \param[in,out] B */
630 /* > \verbatim */
631 /* >          B is COMPLEX array, dimension (LDB, N) */
632 /* >          On entry, the Hermitian matrix B.  If UPLO = 'U', the */
633 /* >          leading N-by-N upper triangular part of B contains the */
634 /* >          upper triangular part of the matrix B.  If UPLO = 'L', */
635 /* >          the leading N-by-N lower triangular part of B contains */
636 /* >          the lower triangular part of the matrix B. */
637 /* > */
638 /* >          On exit, if INFO <= N, the part of B containing the matrix is */
639 /* >          overwritten by the triangular factor U or L from the Cholesky */
640 /* >          factorization B = U**H*U or B = L*L**H. */
641 /* > \endverbatim */
642 /* > */
643 /* > \param[in] LDB */
644 /* > \verbatim */
645 /* >          LDB is INTEGER */
646 /* >          The leading dimension of the array B.  LDB >= f2cmax(1,N). */
647 /* > \endverbatim */
648 /* > */
649 /* > \param[in] VL */
650 /* > \verbatim */
651 /* >          VL is REAL */
652 /* > */
653 /* >          If RANGE='V', the lower bound of the interval to */
654 /* >          be searched for eigenvalues. VL < VU. */
655 /* >          Not referenced if RANGE = 'A' or 'I'. */
656 /* > \endverbatim */
657 /* > */
658 /* > \param[in] VU */
659 /* > \verbatim */
660 /* >          VU is REAL */
661 /* > */
662 /* >          If RANGE='V', the upper bound of the interval to */
663 /* >          be searched for eigenvalues. VL < VU. */
664 /* >          Not referenced if RANGE = 'A' or 'I'. */
665 /* > \endverbatim */
666 /* > */
667 /* > \param[in] IL */
668 /* > \verbatim */
669 /* >          IL is INTEGER */
670 /* > */
671 /* >          If RANGE='I', the index of the */
672 /* >          smallest eigenvalue to be returned. */
673 /* >          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. */
674 /* >          Not referenced if RANGE = 'A' or 'V'. */
675 /* > \endverbatim */
676 /* > */
677 /* > \param[in] IU */
678 /* > \verbatim */
679 /* >          IU is INTEGER */
680 /* > */
681 /* >          If RANGE='I', the index of the */
682 /* >          largest eigenvalue to be returned. */
683 /* >          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. */
684 /* >          Not referenced if RANGE = 'A' or 'V'. */
685 /* > \endverbatim */
686 /* > */
687 /* > \param[in] ABSTOL */
688 /* > \verbatim */
689 /* >          ABSTOL is REAL */
690 /* >          The absolute error tolerance for the eigenvalues. */
691 /* >          An approximate eigenvalue is accepted as converged */
692 /* >          when it is determined to lie in an interval [a,b] */
693 /* >          of width less than or equal to */
694 /* > */
695 /* >                  ABSTOL + EPS *   f2cmax( |a|,|b| ) , */
696 /* > */
697 /* >          where EPS is the machine precision.  If ABSTOL is less than */
698 /* >          or equal to zero, then  EPS*|T|  will be used in its place, */
699 /* >          where |T| is the 1-norm of the tridiagonal matrix obtained */
700 /* >          by reducing C to tridiagonal form, where C is the symmetric */
701 /* >          matrix of the standard symmetric problem to which the */
702 /* >          generalized problem is transformed. */
703 /* > */
704 /* >          Eigenvalues will be computed most accurately when ABSTOL is */
705 /* >          set to twice the underflow threshold 2*SLAMCH('S'), not zero. */
706 /* >          If this routine returns with INFO>0, indicating that some */
707 /* >          eigenvectors did not converge, try setting ABSTOL to */
708 /* >          2*SLAMCH('S'). */
709 /* > \endverbatim */
710 /* > */
711 /* > \param[out] M */
712 /* > \verbatim */
713 /* >          M is INTEGER */
714 /* >          The total number of eigenvalues found.  0 <= M <= N. */
715 /* >          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1. */
716 /* > \endverbatim */
717 /* > */
718 /* > \param[out] W */
719 /* > \verbatim */
720 /* >          W is REAL array, dimension (N) */
721 /* >          The first M elements contain the selected */
722 /* >          eigenvalues in ascending order. */
723 /* > \endverbatim */
724 /* > */
725 /* > \param[out] Z */
726 /* > \verbatim */
727 /* >          Z is COMPLEX array, dimension (LDZ, f2cmax(1,M)) */
728 /* >          If JOBZ = 'N', then Z is not referenced. */
729 /* >          If JOBZ = 'V', then if INFO = 0, the first M columns of Z */
730 /* >          contain the orthonormal eigenvectors of the matrix A */
731 /* >          corresponding to the selected eigenvalues, with the i-th */
732 /* >          column of Z holding the eigenvector associated with W(i). */
733 /* >          The eigenvectors are normalized as follows: */
734 /* >          if ITYPE = 1 or 2, Z**T*B*Z = I; */
735 /* >          if ITYPE = 3, Z**T*inv(B)*Z = I. */
736 /* > */
737 /* >          If an eigenvector fails to converge, then that column of Z */
738 /* >          contains the latest approximation to the eigenvector, and the */
739 /* >          index of the eigenvector is returned in IFAIL. */
740 /* >          Note: the user must ensure that at least f2cmax(1,M) columns are */
741 /* >          supplied in the array Z; if RANGE = 'V', the exact value of M */
742 /* >          is not known in advance and an upper bound must be used. */
743 /* > \endverbatim */
744 /* > */
745 /* > \param[in] LDZ */
746 /* > \verbatim */
747 /* >          LDZ is INTEGER */
748 /* >          The leading dimension of the array Z.  LDZ >= 1, and if */
749 /* >          JOBZ = 'V', LDZ >= f2cmax(1,N). */
750 /* > \endverbatim */
751 /* > */
752 /* > \param[out] WORK */
753 /* > \verbatim */
754 /* >          WORK is COMPLEX array, dimension (MAX(1,LWORK)) */
755 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
756 /* > \endverbatim */
757 /* > */
758 /* > \param[in] LWORK */
759 /* > \verbatim */
760 /* >          LWORK is INTEGER */
761 /* >          The length of the array WORK.  LWORK >= f2cmax(1,2*N). */
762 /* >          For optimal efficiency, LWORK >= (NB+1)*N, */
763 /* >          where NB is the blocksize for CHETRD returned by ILAENV. */
764 /* > */
765 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
766 /* >          only calculates the optimal size of the WORK array, returns */
767 /* >          this value as the first entry of the WORK array, and no error */
768 /* >          message related to LWORK is issued by XERBLA. */
769 /* > \endverbatim */
770 /* > */
771 /* > \param[out] RWORK */
772 /* > \verbatim */
773 /* >          RWORK is REAL array, dimension (7*N) */
774 /* > \endverbatim */
775 /* > */
776 /* > \param[out] IWORK */
777 /* > \verbatim */
778 /* >          IWORK is INTEGER array, dimension (5*N) */
779 /* > \endverbatim */
780 /* > */
781 /* > \param[out] IFAIL */
782 /* > \verbatim */
783 /* >          IFAIL is INTEGER array, dimension (N) */
784 /* >          If JOBZ = 'V', then if INFO = 0, the first M elements of */
785 /* >          IFAIL are zero.  If INFO > 0, then IFAIL contains the */
786 /* >          indices of the eigenvectors that failed to converge. */
787 /* >          If JOBZ = 'N', then IFAIL is not referenced. */
788 /* > \endverbatim */
789 /* > */
790 /* > \param[out] INFO */
791 /* > \verbatim */
792 /* >          INFO is INTEGER */
793 /* >          = 0:  successful exit */
794 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value */
795 /* >          > 0:  CPOTRF or CHEEVX returned an error code: */
796 /* >             <= N:  if INFO = i, CHEEVX failed to converge; */
797 /* >                    i eigenvectors failed to converge.  Their indices */
798 /* >                    are stored in array IFAIL. */
799 /* >             > N:   if INFO = N + i, for 1 <= i <= N, then the leading */
800 /* >                    minor of order i of B is not positive definite. */
801 /* >                    The factorization of B could not be completed and */
802 /* >                    no eigenvalues or eigenvectors were computed. */
803 /* > \endverbatim */
804
805 /*  Authors: */
806 /*  ======== */
807
808 /* > \author Univ. of Tennessee */
809 /* > \author Univ. of California Berkeley */
810 /* > \author Univ. of Colorado Denver */
811 /* > \author NAG Ltd. */
812
813 /* > \date June 2016 */
814
815 /* > \ingroup complexHEeigen */
816
817 /* > \par Contributors: */
818 /*  ================== */
819 /* > */
820 /* >     Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA */
821
822 /*  ===================================================================== */
823 /* Subroutine */ int chegvx_(integer *itype, char *jobz, char *range, char *
824         uplo, integer *n, complex *a, integer *lda, complex *b, integer *ldb, 
825         real *vl, real *vu, integer *il, integer *iu, real *abstol, integer *
826         m, real *w, complex *z__, integer *ldz, complex *work, integer *lwork,
827          real *rwork, integer *iwork, integer *ifail, integer *info)
828 {
829     /* System generated locals */
830     integer a_dim1, a_offset, b_dim1, b_offset, z_dim1, z_offset, i__1, i__2;
831
832     /* Local variables */
833     extern logical lsame_(char *, char *);
834     extern /* Subroutine */ int ctrmm_(char *, char *, char *, char *, 
835             integer *, integer *, complex *, complex *, integer *, complex *, 
836             integer *);
837     char trans[1];
838     extern /* Subroutine */ int ctrsm_(char *, char *, char *, char *, 
839             integer *, integer *, complex *, complex *, integer *, complex *, 
840             integer *);
841     logical upper, wantz;
842     integer nb;
843     logical alleig, indeig, valeig;
844     extern /* Subroutine */ int chegst_(integer *, char *, integer *, complex 
845             *, integer *, complex *, integer *, integer *);
846     extern integer ilaenv_(integer *, char *, char *, integer *, integer *, 
847             integer *, integer *, ftnlen, ftnlen);
848     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen), cheevx_(
849             char *, char *, char *, integer *, complex *, integer *, real *, 
850             real *, integer *, integer *, real *, integer *, real *, complex *
851             , integer *, complex *, integer *, real *, integer *, integer *, 
852             integer *), cpotrf_(char *, integer *, 
853             complex *, integer *, integer *);
854     integer lwkopt;
855     logical lquery;
856
857
858 /*  -- LAPACK driver routine (version 3.7.0) -- */
859 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
860 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
861 /*     June 2016 */
862
863
864 /*  ===================================================================== */
865
866
867 /*     Test the input parameters. */
868
869     /* Parameter adjustments */
870     a_dim1 = *lda;
871     a_offset = 1 + a_dim1 * 1;
872     a -= a_offset;
873     b_dim1 = *ldb;
874     b_offset = 1 + b_dim1 * 1;
875     b -= b_offset;
876     --w;
877     z_dim1 = *ldz;
878     z_offset = 1 + z_dim1 * 1;
879     z__ -= z_offset;
880     --work;
881     --rwork;
882     --iwork;
883     --ifail;
884
885     /* Function Body */
886     wantz = lsame_(jobz, "V");
887     upper = lsame_(uplo, "U");
888     alleig = lsame_(range, "A");
889     valeig = lsame_(range, "V");
890     indeig = lsame_(range, "I");
891     lquery = *lwork == -1;
892
893     *info = 0;
894     if (*itype < 1 || *itype > 3) {
895         *info = -1;
896     } else if (! (wantz || lsame_(jobz, "N"))) {
897         *info = -2;
898     } else if (! (alleig || valeig || indeig)) {
899         *info = -3;
900     } else if (! (upper || lsame_(uplo, "L"))) {
901         *info = -4;
902     } else if (*n < 0) {
903         *info = -5;
904     } else if (*lda < f2cmax(1,*n)) {
905         *info = -7;
906     } else if (*ldb < f2cmax(1,*n)) {
907         *info = -9;
908     } else {
909         if (valeig) {
910             if (*n > 0 && *vu <= *vl) {
911                 *info = -11;
912             }
913         } else if (indeig) {
914             if (*il < 1 || *il > f2cmax(1,*n)) {
915                 *info = -12;
916             } else if (*iu < f2cmin(*n,*il) || *iu > *n) {
917                 *info = -13;
918             }
919         }
920     }
921     if (*info == 0) {
922         if (*ldz < 1 || wantz && *ldz < *n) {
923             *info = -18;
924         }
925     }
926
927     if (*info == 0) {
928         nb = ilaenv_(&c__1, "CHETRD", uplo, n, &c_n1, &c_n1, &c_n1, (ftnlen)6,
929                  (ftnlen)1);
930 /* Computing MAX */
931         i__1 = 1, i__2 = (nb + 1) * *n;
932         lwkopt = f2cmax(i__1,i__2);
933         work[1].r = (real) lwkopt, work[1].i = 0.f;
934
935 /* Computing MAX */
936         i__1 = 1, i__2 = *n << 1;
937         if (*lwork < f2cmax(i__1,i__2) && ! lquery) {
938             *info = -20;
939         }
940     }
941
942     if (*info != 0) {
943         i__1 = -(*info);
944         xerbla_("CHEGVX", &i__1, (ftnlen)6);
945         return 0;
946     } else if (lquery) {
947         return 0;
948     }
949
950 /*     Quick return if possible */
951
952     *m = 0;
953     if (*n == 0) {
954         return 0;
955     }
956
957 /*     Form a Cholesky factorization of B. */
958
959     cpotrf_(uplo, n, &b[b_offset], ldb, info);
960     if (*info != 0) {
961         *info = *n + *info;
962         return 0;
963     }
964
965 /*     Transform problem to standard eigenvalue problem and solve. */
966
967     chegst_(itype, uplo, n, &a[a_offset], lda, &b[b_offset], ldb, info);
968     cheevx_(jobz, range, uplo, n, &a[a_offset], lda, vl, vu, il, iu, abstol, 
969             m, &w[1], &z__[z_offset], ldz, &work[1], lwork, &rwork[1], &iwork[
970             1], &ifail[1], info);
971
972     if (wantz) {
973
974 /*        Backtransform eigenvectors to the original problem. */
975
976         if (*info > 0) {
977             *m = *info - 1;
978         }
979         if (*itype == 1 || *itype == 2) {
980
981 /*           For A*x=(lambda)*B*x and A*B*x=(lambda)*x; */
982 /*           backtransform eigenvectors: x = inv(L)**H*y or inv(U)*y */
983
984             if (upper) {
985                 *(unsigned char *)trans = 'N';
986             } else {
987                 *(unsigned char *)trans = 'C';
988             }
989
990             ctrsm_("Left", uplo, trans, "Non-unit", n, m, &c_b1, &b[b_offset],
991                      ldb, &z__[z_offset], ldz);
992
993         } else if (*itype == 3) {
994
995 /*           For B*A*x=(lambda)*x; */
996 /*           backtransform eigenvectors: x = L*y or U**H*y */
997
998             if (upper) {
999                 *(unsigned char *)trans = 'C';
1000             } else {
1001                 *(unsigned char *)trans = 'N';
1002             }
1003
1004             ctrmm_("Left", uplo, trans, "Non-unit", n, m, &c_b1, &b[b_offset],
1005                      ldb, &z__[z_offset], ldz);
1006         }
1007     }
1008
1009 /*     Set WORK(1) to optimal complex workspace size. */
1010
1011     work[1].r = (real) lwkopt, work[1].i = 0.f;
1012
1013     return 0;
1014
1015 /*     End of CHEGVX */
1016
1017 } /* chegvx_ */
1018