C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / chbev_2stage.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__2 = 2;
516 static integer c_n1 = -1;
517 static integer c__3 = 3;
518 static integer c__4 = 4;
519 static real c_b21 = 1.f;
520 static integer c__1 = 1;
521
522 /* > \brief <b> CHBEV_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for 
523 OTHER matrices</b> */
524
525 /*  @generated from zhbev_2stage.f, fortran z -> c, Sat Nov  5 23:18:20 2016 */
526
527 /*  =========== DOCUMENTATION =========== */
528
529 /* Online html documentation available at */
530 /*            http://www.netlib.org/lapack/explore-html/ */
531
532 /* > \htmlonly */
533 /* > Download CHBEV_2STAGE + dependencies */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/chbev_2
535 stage.f"> */
536 /* > [TGZ]</a> */
537 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/chbev_2
538 stage.f"> */
539 /* > [ZIP]</a> */
540 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/chbev_2
541 stage.f"> */
542 /* > [TXT]</a> */
543 /* > \endhtmlonly */
544
545 /*  Definition: */
546 /*  =========== */
547
548 /*       SUBROUTINE CHBEV_2STAGE( JOBZ, UPLO, N, KD, AB, LDAB, W, Z, LDZ, */
549 /*                                WORK, LWORK, RWORK, INFO ) */
550
551 /*       IMPLICIT NONE */
552
553 /*       CHARACTER          JOBZ, UPLO */
554 /*       INTEGER            INFO, KD, LDAB, LDZ, N, LWORK */
555 /*       REAL               RWORK( * ), W( * ) */
556 /*       COMPLEX            AB( LDAB, * ), WORK( * ), Z( LDZ, * ) */
557
558
559 /* > \par Purpose: */
560 /*  ============= */
561 /* > */
562 /* > \verbatim */
563 /* > */
564 /* > CHBEV_2STAGE computes all the eigenvalues and, optionally, eigenvectors of */
565 /* > a complex Hermitian band matrix A using the 2stage technique for */
566 /* > the reduction to tridiagonal. */
567 /* > \endverbatim */
568
569 /*  Arguments: */
570 /*  ========== */
571
572 /* > \param[in] JOBZ */
573 /* > \verbatim */
574 /* >          JOBZ is CHARACTER*1 */
575 /* >          = 'N':  Compute eigenvalues only; */
576 /* >          = 'V':  Compute eigenvalues and eigenvectors. */
577 /* >                  Not available in this release. */
578 /* > \endverbatim */
579 /* > */
580 /* > \param[in] UPLO */
581 /* > \verbatim */
582 /* >          UPLO is CHARACTER*1 */
583 /* >          = 'U':  Upper triangle of A is stored; */
584 /* >          = 'L':  Lower triangle of A is stored. */
585 /* > \endverbatim */
586 /* > */
587 /* > \param[in] N */
588 /* > \verbatim */
589 /* >          N is INTEGER */
590 /* >          The order of the matrix A.  N >= 0. */
591 /* > \endverbatim */
592 /* > */
593 /* > \param[in] KD */
594 /* > \verbatim */
595 /* >          KD is INTEGER */
596 /* >          The number of superdiagonals of the matrix A if UPLO = 'U', */
597 /* >          or the number of subdiagonals if UPLO = 'L'.  KD >= 0. */
598 /* > \endverbatim */
599 /* > */
600 /* > \param[in,out] AB */
601 /* > \verbatim */
602 /* >          AB is COMPLEX array, dimension (LDAB, N) */
603 /* >          On entry, the upper or lower triangle of the Hermitian band */
604 /* >          matrix A, stored in the first KD+1 rows of the array.  The */
605 /* >          j-th column of A is stored in the j-th column of the array AB */
606 /* >          as follows: */
607 /* >          if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for f2cmax(1,j-kd)<=i<=j; */
608 /* >          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=f2cmin(n,j+kd). */
609 /* > */
610 /* >          On exit, AB is overwritten by values generated during the */
611 /* >          reduction to tridiagonal form.  If UPLO = 'U', the first */
612 /* >          superdiagonal and the diagonal of the tridiagonal matrix T */
613 /* >          are returned in rows KD and KD+1 of AB, and if UPLO = 'L', */
614 /* >          the diagonal and first subdiagonal of T are returned in the */
615 /* >          first two rows of AB. */
616 /* > \endverbatim */
617 /* > */
618 /* > \param[in] LDAB */
619 /* > \verbatim */
620 /* >          LDAB is INTEGER */
621 /* >          The leading dimension of the array AB.  LDAB >= KD + 1. */
622 /* > \endverbatim */
623 /* > */
624 /* > \param[out] W */
625 /* > \verbatim */
626 /* >          W is REAL array, dimension (N) */
627 /* >          If INFO = 0, the eigenvalues in ascending order. */
628 /* > \endverbatim */
629 /* > */
630 /* > \param[out] Z */
631 /* > \verbatim */
632 /* >          Z is COMPLEX array, dimension (LDZ, N) */
633 /* >          If JOBZ = 'V', then if INFO = 0, Z contains the orthonormal */
634 /* >          eigenvectors of the matrix A, with the i-th column of Z */
635 /* >          holding the eigenvector associated with W(i). */
636 /* >          If JOBZ = 'N', then Z is not referenced. */
637 /* > \endverbatim */
638 /* > */
639 /* > \param[in] LDZ */
640 /* > \verbatim */
641 /* >          LDZ is INTEGER */
642 /* >          The leading dimension of the array Z.  LDZ >= 1, and if */
643 /* >          JOBZ = 'V', LDZ >= f2cmax(1,N). */
644 /* > \endverbatim */
645 /* > */
646 /* > \param[out] WORK */
647 /* > \verbatim */
648 /* >          WORK is COMPLEX array, dimension LWORK */
649 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
650 /* > \endverbatim */
651 /* > */
652 /* > \param[in] LWORK */
653 /* > \verbatim */
654 /* >          LWORK is INTEGER */
655 /* >          The length of the array WORK. LWORK >= 1, when N <= 1; */
656 /* >          otherwise */
657 /* >          If JOBZ = 'N' and N > 1, LWORK must be queried. */
658 /* >                                   LWORK = MAX(1, dimension) where */
659 /* >                                   dimension = (2KD+1)*N + KD*NTHREADS */
660 /* >                                   where KD is the size of the band. */
661 /* >                                   NTHREADS is the number of threads used when */
662 /* >                                   openMP compilation is enabled, otherwise =1. */
663 /* >          If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available. */
664 /* > */
665 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
666 /* >          only calculates the optimal sizes of the WORK, RWORK and */
667 /* >          IWORK arrays, returns these values as the first entries of */
668 /* >          the WORK, RWORK and IWORK arrays, and no error message */
669 /* >          related to LWORK or LRWORK or LIWORK is issued by XERBLA. */
670 /* > \endverbatim */
671 /* > */
672 /* > \param[out] RWORK */
673 /* > \verbatim */
674 /* >          RWORK is REAL array, dimension (f2cmax(1,3*N-2)) */
675 /* > \endverbatim */
676 /* > */
677 /* > \param[out] INFO */
678 /* > \verbatim */
679 /* >          INFO is INTEGER */
680 /* >          = 0:  successful exit. */
681 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value. */
682 /* >          > 0:  if INFO = i, the algorithm failed to converge; i */
683 /* >                off-diagonal elements of an intermediate tridiagonal */
684 /* >                form did not converge to zero. */
685 /* > \endverbatim */
686
687 /*  Authors: */
688 /*  ======== */
689
690 /* > \author Univ. of Tennessee */
691 /* > \author Univ. of California Berkeley */
692 /* > \author Univ. of Colorado Denver */
693 /* > \author NAG Ltd. */
694
695 /* > \date November 2017 */
696
697 /* > \ingroup complexOTHEReigen */
698
699 /* > \par Further Details: */
700 /*  ===================== */
701 /* > */
702 /* > \verbatim */
703 /* > */
704 /* >  All details about the 2stage techniques are available in: */
705 /* > */
706 /* >  Azzam Haidar, Hatem Ltaief, and Jack Dongarra. */
707 /* >  Parallel reduction to condensed forms for symmetric eigenvalue problems */
708 /* >  using aggregated fine-grained and memory-aware kernels. In Proceedings */
709 /* >  of 2011 International Conference for High Performance Computing, */
710 /* >  Networking, Storage and Analysis (SC '11), New York, NY, USA, */
711 /* >  Article 8 , 11 pages. */
712 /* >  http://doi.acm.org/10.1145/2063384.2063394 */
713 /* > */
714 /* >  A. Haidar, J. Kurzak, P. Luszczek, 2013. */
715 /* >  An improved parallel singular value algorithm and its implementation */
716 /* >  for multicore hardware, In Proceedings of 2013 International Conference */
717 /* >  for High Performance Computing, Networking, Storage and Analysis (SC '13). */
718 /* >  Denver, Colorado, USA, 2013. */
719 /* >  Article 90, 12 pages. */
720 /* >  http://doi.acm.org/10.1145/2503210.2503292 */
721 /* > */
722 /* >  A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra. */
723 /* >  A novel hybrid CPU-GPU generalized eigensolver for electronic structure */
724 /* >  calculations based on fine-grained memory aware tasks. */
725 /* >  International Journal of High Performance Computing Applications. */
726 /* >  Volume 28 Issue 2, Pages 196-209, May 2014. */
727 /* >  http://hpc.sagepub.com/content/28/2/196 */
728 /* > */
729 /* > \endverbatim */
730
731 /*  ===================================================================== */
732 /* Subroutine */ int chbev_2stage_(char *jobz, char *uplo, integer *n, 
733         integer *kd, complex *ab, integer *ldab, real *w, complex *z__, 
734         integer *ldz, complex *work, integer *lwork, real *rwork, integer *
735         info)
736 {
737     /* System generated locals */
738     integer ab_dim1, ab_offset, z_dim1, z_offset, i__1;
739     real r__1;
740
741     /* Local variables */
742     extern /* Subroutine */ int chetrd_hb2st_(char *, char *, char *, 
743             integer *, integer *, complex *, integer *, real *, real *, 
744             complex *, integer *, complex *, integer *, integer *);
745     integer inde;
746     extern integer ilaenv2stage_(integer *, char *, char *, integer *, 
747             integer *, integer *, integer *);
748     real anrm;
749     integer imax;
750     real rmin, rmax, sigma;
751     extern logical lsame_(char *, char *);
752     integer iinfo;
753     extern /* Subroutine */ int sscal_(integer *, real *, real *, integer *);
754     integer lhtrd, lwmin;
755     logical lower;
756     integer lwtrd;
757     logical wantz;
758     integer ib;
759     extern real clanhb_(char *, char *, integer *, integer *, complex *, 
760             integer *, real *);
761     integer iscale;
762     extern /* Subroutine */ int clascl_(char *, integer *, integer *, real *, 
763             real *, integer *, integer *, complex *, integer *, integer *);
764     extern real slamch_(char *);
765     real safmin;
766     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
767     real bignum;
768     integer indwrk, indrwk;
769     extern /* Subroutine */ int csteqr_(char *, integer *, real *, real *, 
770             complex *, integer *, real *, integer *), ssterf_(integer 
771             *, real *, real *, integer *);
772     integer llwork;
773     real smlnum;
774     logical lquery;
775     real eps;
776     integer indhous;
777
778
779
780 /*  -- LAPACK driver routine (version 3.8.0) -- */
781 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
782 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
783 /*     November 2017 */
784
785
786 /*  ===================================================================== */
787
788
789 /*     Test the input parameters. */
790
791     /* Parameter adjustments */
792     ab_dim1 = *ldab;
793     ab_offset = 1 + ab_dim1 * 1;
794     ab -= ab_offset;
795     --w;
796     z_dim1 = *ldz;
797     z_offset = 1 + z_dim1 * 1;
798     z__ -= z_offset;
799     --work;
800     --rwork;
801
802     /* Function Body */
803     wantz = lsame_(jobz, "V");
804     lower = lsame_(uplo, "L");
805     lquery = *lwork == -1;
806
807     *info = 0;
808     if (! lsame_(jobz, "N")) {
809         *info = -1;
810     } else if (! (lower || lsame_(uplo, "U"))) {
811         *info = -2;
812     } else if (*n < 0) {
813         *info = -3;
814     } else if (*kd < 0) {
815         *info = -4;
816     } else if (*ldab < *kd + 1) {
817         *info = -6;
818     } else if (*ldz < 1 || wantz && *ldz < *n) {
819         *info = -9;
820     }
821
822     if (*info == 0) {
823         if (*n <= 1) {
824             lwmin = 1;
825             work[1].r = (real) lwmin, work[1].i = 0.f;
826         } else {
827             ib = ilaenv2stage_(&c__2, "CHETRD_HB2ST", jobz, n, kd, &c_n1, &
828                     c_n1);
829             lhtrd = ilaenv2stage_(&c__3, "CHETRD_HB2ST", jobz, n, kd, &ib, &
830                     c_n1);
831             lwtrd = ilaenv2stage_(&c__4, "CHETRD_HB2ST", jobz, n, kd, &ib, &
832                     c_n1);
833             lwmin = lhtrd + lwtrd;
834             work[1].r = (real) lwmin, work[1].i = 0.f;
835         }
836
837         if (*lwork < lwmin && ! lquery) {
838             *info = -11;
839         }
840     }
841
842     if (*info != 0) {
843         i__1 = -(*info);
844         xerbla_("CHBEV_2STAGE ", &i__1, (ftnlen)13);
845         return 0;
846     } else if (lquery) {
847         return 0;
848     }
849
850 /*     Quick return if possible */
851
852     if (*n == 0) {
853         return 0;
854     }
855
856     if (*n == 1) {
857         if (lower) {
858             i__1 = ab_dim1 + 1;
859             w[1] = ab[i__1].r;
860         } else {
861             i__1 = *kd + 1 + ab_dim1;
862             w[1] = ab[i__1].r;
863         }
864         if (wantz) {
865             i__1 = z_dim1 + 1;
866             z__[i__1].r = 1.f, z__[i__1].i = 0.f;
867         }
868         return 0;
869     }
870
871 /*     Get machine constants. */
872
873     safmin = slamch_("Safe minimum");
874     eps = slamch_("Precision");
875     smlnum = safmin / eps;
876     bignum = 1.f / smlnum;
877     rmin = sqrt(smlnum);
878     rmax = sqrt(bignum);
879
880 /*     Scale matrix to allowable range, if necessary. */
881
882     anrm = clanhb_("M", uplo, n, kd, &ab[ab_offset], ldab, &rwork[1]);
883     iscale = 0;
884     if (anrm > 0.f && anrm < rmin) {
885         iscale = 1;
886         sigma = rmin / anrm;
887     } else if (anrm > rmax) {
888         iscale = 1;
889         sigma = rmax / anrm;
890     }
891     if (iscale == 1) {
892         if (lower) {
893             clascl_("B", kd, kd, &c_b21, &sigma, n, n, &ab[ab_offset], ldab, 
894                     info);
895         } else {
896             clascl_("Q", kd, kd, &c_b21, &sigma, n, n, &ab[ab_offset], ldab, 
897                     info);
898         }
899     }
900
901 /*     Call CHBTRD_HB2ST to reduce Hermitian band matrix to tridiagonal form. */
902
903     inde = 1;
904     indhous = 1;
905     indwrk = indhous + lhtrd;
906     llwork = *lwork - indwrk + 1;
907
908     chetrd_hb2st_("N", jobz, uplo, n, kd, &ab[ab_offset], ldab, &w[1], &
909             rwork[inde], &work[indhous], &lhtrd, &work[indwrk], &llwork, &
910             iinfo);
911
912 /*     For eigenvalues only, call SSTERF.  For eigenvectors, call CSTEQR. */
913
914     if (! wantz) {
915         ssterf_(n, &w[1], &rwork[inde], info);
916     } else {
917         indrwk = inde + *n;
918         csteqr_(jobz, n, &w[1], &rwork[inde], &z__[z_offset], ldz, &rwork[
919                 indrwk], info);
920     }
921
922 /*     If matrix was scaled, then rescale eigenvalues appropriately. */
923
924     if (iscale == 1) {
925         if (*info == 0) {
926             imax = *n;
927         } else {
928             imax = *info - 1;
929         }
930         r__1 = 1.f / sigma;
931         sscal_(&imax, &r__1, &w[1], &c__1);
932     }
933
934 /*     Set WORK(1) to optimal workspace size. */
935
936     work[1].r = (real) lwmin, work[1].i = 0.f;
937
938     return 0;
939
940 /*     End of CHBEV_2STAGE */
941
942 } /* chbev_2stage__ */
943