C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / chetrd_hb2st.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 = {0.f,0.f};
516 static integer c__2 = 2;
517 static integer c_n1 = -1;
518 static integer c__3 = 3;
519 static integer c__4 = 4;
520
521 /* > \brief \b CHBTRD_HB2ST reduces a complex Hermitian band matrix A to real symmetric tridiagonal form T */
522
523 /*  =========== DOCUMENTATION =========== */
524
525 /* Online html documentation available at */
526 /*            http://www.netlib.org/lapack/explore-html/ */
527
528 /* > \htmlonly */
529 /* > Download CHBTRD_HB2ST + dependencies */
530 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/chbtrd_
531 hb2st.f"> */
532 /* > [TGZ]</a> */
533 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/chbtrd_
534 hb2st.f"> */
535 /* > [ZIP]</a> */
536 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/chbtrd_
537 hb2st.f"> */
538 /* > [TXT]</a> */
539 /* > \endhtmlonly */
540
541 /*  Definition: */
542 /*  =========== */
543
544 /*       SUBROUTINE CHETRD_HB2ST( STAGE1, VECT, UPLO, N, KD, AB, LDAB, */
545 /*                               D, E, HOUS, LHOUS, WORK, LWORK, INFO ) */
546
547 /*       #if defined(_OPENMP) */
548 /*       use omp_lib */
549 /*       #endif */
550
551 /*       IMPLICIT NONE */
552
553 /*       CHARACTER          STAGE1, UPLO, VECT */
554 /*       INTEGER            N, KD, IB, LDAB, LHOUS, LWORK, INFO */
555 /*       REAL               D( * ), E( * ) */
556 /*       COMPLEX            AB( LDAB, * ), HOUS( * ), WORK( * ) */
557
558
559 /* > \par Purpose: */
560 /*  ============= */
561 /* > */
562 /* > \verbatim */
563 /* > */
564 /* > CHETRD_HB2ST reduces a complex Hermitian band matrix A to real symmetric */
565 /* > tridiagonal form T by a unitary similarity transformation: */
566 /* > Q**H * A * Q = T. */
567 /* > \endverbatim */
568
569 /*  Arguments: */
570 /*  ========== */
571
572 /* > \param[in] STAGE1 */
573 /* > \verbatim */
574 /* >          STAGE1 is CHARACTER*1 */
575 /* >          = 'N':  "No": to mention that the stage 1 of the reduction */
576 /* >                  from dense to band using the chetrd_he2hb routine */
577 /* >                  was not called before this routine to reproduce AB. */
578 /* >                  In other term this routine is called as standalone. */
579 /* >          = 'Y':  "Yes": to mention that the stage 1 of the */
580 /* >                  reduction from dense to band using the chetrd_he2hb */
581 /* >                  routine has been called to produce AB (e.g., AB is */
582 /* >                  the output of chetrd_he2hb. */
583 /* > \endverbatim */
584 /* > */
585 /* > \param[in] VECT */
586 /* > \verbatim */
587 /* >          VECT is CHARACTER*1 */
588 /* >          = 'N':  No need for the Housholder representation, */
589 /* >                  and thus LHOUS is of size f2cmax(1, 4*N); */
590 /* >          = 'V':  the Householder representation is needed to */
591 /* >                  either generate or to apply Q later on, */
592 /* >                  then LHOUS is to be queried and computed. */
593 /* >                  (NOT AVAILABLE IN THIS RELEASE). */
594 /* > \endverbatim */
595 /* > */
596 /* > \param[in] UPLO */
597 /* > \verbatim */
598 /* >          UPLO is CHARACTER*1 */
599 /* >          = 'U':  Upper triangle of A is stored; */
600 /* >          = 'L':  Lower triangle of A is stored. */
601 /* > \endverbatim */
602 /* > */
603 /* > \param[in] N */
604 /* > \verbatim */
605 /* >          N is INTEGER */
606 /* >          The order of the matrix A.  N >= 0. */
607 /* > \endverbatim */
608 /* > */
609 /* > \param[in] KD */
610 /* > \verbatim */
611 /* >          KD is INTEGER */
612 /* >          The number of superdiagonals of the matrix A if UPLO = 'U', */
613 /* >          or the number of subdiagonals if UPLO = 'L'.  KD >= 0. */
614 /* > \endverbatim */
615 /* > */
616 /* > \param[in,out] AB */
617 /* > \verbatim */
618 /* >          AB is COMPLEX array, dimension (LDAB,N) */
619 /* >          On entry, the upper or lower triangle of the Hermitian band */
620 /* >          matrix A, stored in the first KD+1 rows of the array.  The */
621 /* >          j-th column of A is stored in the j-th column of the array AB */
622 /* >          as follows: */
623 /* >          if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for f2cmax(1,j-kd)<=i<=j; */
624 /* >          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=f2cmin(n,j+kd). */
625 /* >          On exit, the diagonal elements of AB are overwritten by the */
626 /* >          diagonal elements of the tridiagonal matrix T; if KD > 0, the */
627 /* >          elements on the first superdiagonal (if UPLO = 'U') or the */
628 /* >          first subdiagonal (if UPLO = 'L') are overwritten by the */
629 /* >          off-diagonal elements of T; the rest of AB is overwritten by */
630 /* >          values generated during the reduction. */
631 /* > \endverbatim */
632 /* > */
633 /* > \param[in] LDAB */
634 /* > \verbatim */
635 /* >          LDAB is INTEGER */
636 /* >          The leading dimension of the array AB.  LDAB >= KD+1. */
637 /* > \endverbatim */
638 /* > */
639 /* > \param[out] D */
640 /* > \verbatim */
641 /* >          D is REAL array, dimension (N) */
642 /* >          The diagonal elements of the tridiagonal matrix T. */
643 /* > \endverbatim */
644 /* > */
645 /* > \param[out] E */
646 /* > \verbatim */
647 /* >          E is REAL array, dimension (N-1) */
648 /* >          The off-diagonal elements of the tridiagonal matrix T: */
649 /* >          E(i) = T(i,i+1) if UPLO = 'U'; E(i) = T(i+1,i) if UPLO = 'L'. */
650 /* > \endverbatim */
651 /* > */
652 /* > \param[out] HOUS */
653 /* > \verbatim */
654 /* >          HOUS is COMPLEX array, dimension LHOUS, that */
655 /* >          store the Householder representation. */
656 /* > \endverbatim */
657 /* > */
658 /* > \param[in] LHOUS */
659 /* > \verbatim */
660 /* >          LHOUS is INTEGER */
661 /* >          The dimension of the array HOUS. LHOUS = MAX(1, dimension) */
662 /* >          If LWORK = -1, or LHOUS=-1, */
663 /* >          then a query is assumed; the routine */
664 /* >          only calculates the optimal size of the HOUS array, returns */
665 /* >          this value as the first entry of the HOUS array, and no error */
666 /* >          message related to LHOUS is issued by XERBLA. */
667 /* >          LHOUS = MAX(1, dimension) where */
668 /* >          dimension = 4*N if VECT='N' */
669 /* >          not available now if VECT='H' */
670 /* > \endverbatim */
671 /* > */
672 /* > \param[out] WORK */
673 /* > \verbatim */
674 /* >          WORK is COMPLEX array, dimension LWORK. */
675 /* > \endverbatim */
676 /* > */
677 /* > \param[in] LWORK */
678 /* > \verbatim */
679 /* >          LWORK is INTEGER */
680 /* >          The dimension of the array WORK. LWORK = MAX(1, dimension) */
681 /* >          If LWORK = -1, or LHOUS=-1, */
682 /* >          then a workspace query is assumed; the routine */
683 /* >          only calculates the optimal size of the WORK array, returns */
684 /* >          this value as the first entry of the WORK array, and no error */
685 /* >          message related to LWORK is issued by XERBLA. */
686 /* >          LWORK = MAX(1, dimension) where */
687 /* >          dimension   = (2KD+1)*N + KD*NTHREADS */
688 /* >          where KD is the blocking size of the reduction, */
689 /* >          FACTOPTNB is the blocking used by the QR or LQ */
690 /* >          algorithm, usually FACTOPTNB=128 is a good choice */
691 /* >          NTHREADS is the number of threads used when */
692 /* >          openMP compilation is enabled, otherwise =1. */
693 /* > \endverbatim */
694 /* > */
695 /* > \param[out] INFO */
696 /* > \verbatim */
697 /* >          INFO is INTEGER */
698 /* >          = 0:  successful exit */
699 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value */
700 /* > \endverbatim */
701
702 /*  Authors: */
703 /*  ======== */
704
705 /* > \author Univ. of Tennessee */
706 /* > \author Univ. of California Berkeley */
707 /* > \author Univ. of Colorado Denver */
708 /* > \author NAG Ltd. */
709
710 /* > \date November 2017 */
711
712 /* > \ingroup complexOTHERcomputational */
713
714 /* > \par Further Details: */
715 /*  ===================== */
716 /* > */
717 /* > \verbatim */
718 /* > */
719 /* >  Implemented by Azzam Haidar. */
720 /* > */
721 /* >  All details are available on technical report, SC11, SC13 papers. */
722 /* > */
723 /* >  Azzam Haidar, Hatem Ltaief, and Jack Dongarra. */
724 /* >  Parallel reduction to condensed forms for symmetric eigenvalue problems */
725 /* >  using aggregated fine-grained and memory-aware kernels. In Proceedings */
726 /* >  of 2011 International Conference for High Performance Computing, */
727 /* >  Networking, Storage and Analysis (SC '11), New York, NY, USA, */
728 /* >  Article 8 , 11 pages. */
729 /* >  http://doi.acm.org/10.1145/2063384.2063394 */
730 /* > */
731 /* >  A. Haidar, J. Kurzak, P. Luszczek, 2013. */
732 /* >  An improved parallel singular value algorithm and its implementation */
733 /* >  for multicore hardware, In Proceedings of 2013 International Conference */
734 /* >  for High Performance Computing, Networking, Storage and Analysis (SC '13). */
735 /* >  Denver, Colorado, USA, 2013. */
736 /* >  Article 90, 12 pages. */
737 /* >  http://doi.acm.org/10.1145/2503210.2503292 */
738 /* > */
739 /* >  A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra. */
740 /* >  A novel hybrid CPU-GPU generalized eigensolver for electronic structure */
741 /* >  calculations based on fine-grained memory aware tasks. */
742 /* >  International Journal of High Performance Computing Applications. */
743 /* >  Volume 28 Issue 2, Pages 196-209, May 2014. */
744 /* >  http://hpc.sagepub.com/content/28/2/196 */
745 /* > */
746 /* > \endverbatim */
747 /* > */
748 /*  ===================================================================== */
749 /* Subroutine */ int chetrd_hb2st_(char *stage1, char *vect, char *uplo, 
750         integer *n, integer *kd, complex *ab, integer *ldab, real *d__, real *
751         e, complex *hous, integer *lhous, complex *work, integer *lwork, 
752         integer *info)
753 {
754     /* System generated locals */
755     integer ab_dim1, ab_offset, i__1, i__2, i__3, i__4, i__5;
756     real r__1;
757     complex q__1;
758
759     /* Local variables */
760     integer inda;
761     extern integer ilaenv2stage_(integer *, char *, char *, integer *, 
762             integer *, integer *, integer *);
763     integer thed, indv, myid, indw, apos, dpos, abofdpos, nthreads, i__, k, m,
764              edind, debug;
765     extern logical lsame_(char *, char *);
766     integer lhmin, sicev, sizea, shift, stind, colpt, lwmin, awpos;
767     logical wantq, upper;
768     integer grsiz, ttype, stepercol, ed, ib;
769     extern /* Subroutine */ int chb2st_kernels_(char *, logical *, integer *,
770              integer *, integer *, integer *, integer *, integer *, integer *,
771              complex *, integer *, complex *, complex *, integer *, complex *);
772     integer st, abdpos;
773     extern /* Subroutine */ int clacpy_(char *, integer *, integer *, complex 
774             *, integer *, complex *, integer *), claset_(char *, 
775             integer *, integer *, complex *, complex *, complex *, integer *), xerbla_(char *, integer *, ftnlen);
776     integer thgrid, thgrnb, indtau;
777     real abstmp;
778     integer ofdpos, blklastind;
779     extern /* Subroutine */ int mecago_();
780     logical lquery, afters1;
781     integer lda, tid, ldv;
782     complex tmp;
783     integer stt, sweepid, nbtiles, sizetau, thgrsiz;
784
785
786
787
788
789 /*  -- LAPACK computational routine (version 3.8.0) -- */
790 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
791 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
792 /*     November 2017 */
793
794
795 /*  ===================================================================== */
796
797
798 /*     Determine the minimal workspace size required. */
799 /*     Test the input parameters */
800
801     /* Parameter adjustments */
802     ab_dim1 = *ldab;
803     ab_offset = 1 + ab_dim1 * 1;
804     ab -= ab_offset;
805     --d__;
806     --e;
807     --hous;
808     --work;
809
810     /* Function Body */
811     debug = 0;
812     *info = 0;
813     afters1 = lsame_(stage1, "Y");
814     wantq = lsame_(vect, "V");
815     upper = lsame_(uplo, "U");
816     lquery = *lwork == -1 || *lhous == -1;
817
818 /*     Determine the block size, the workspace size and the hous size. */
819
820     ib = ilaenv2stage_(&c__2, "CHETRD_HB2ST", vect, n, kd, &c_n1, &c_n1);
821     lhmin = ilaenv2stage_(&c__3, "CHETRD_HB2ST", vect, n, kd, &ib, &c_n1);
822     lwmin = ilaenv2stage_(&c__4, "CHETRD_HB2ST", vect, n, kd, &ib, &c_n1);
823
824     if (! afters1 && ! lsame_(stage1, "N")) {
825         *info = -1;
826     } else if (! lsame_(vect, "N")) {
827         *info = -2;
828     } else if (! upper && ! lsame_(uplo, "L")) {
829         *info = -3;
830     } else if (*n < 0) {
831         *info = -4;
832     } else if (*kd < 0) {
833         *info = -5;
834     } else if (*ldab < *kd + 1) {
835         *info = -7;
836     } else if (*lhous < lhmin && ! lquery) {
837         *info = -11;
838     } else if (*lwork < lwmin && ! lquery) {
839         *info = -13;
840     }
841
842     if (*info == 0) {
843         hous[1].r = (real) lhmin, hous[1].i = 0.f;
844         work[1].r = (real) lwmin, work[1].i = 0.f;
845     }
846
847     if (*info != 0) {
848         i__1 = -(*info);
849         xerbla_("CHETRD_HB2ST", &i__1, (ftnlen)12);
850         return 0;
851     } else if (lquery) {
852         return 0;
853     }
854
855 /*     Quick return if possible */
856
857     if (*n == 0) {
858         hous[1].r = 1.f, hous[1].i = 0.f;
859         work[1].r = 1.f, work[1].i = 0.f;
860         return 0;
861     }
862
863 /*     Determine pointer position */
864
865     ldv = *kd + ib;
866     sizetau = *n << 1;
867     sicev = *n << 1;
868     indtau = 1;
869     indv = indtau + sizetau;
870     lda = (*kd << 1) + 1;
871     sizea = lda * *n;
872     inda = 1;
873     indw = inda + sizea;
874     nthreads = 1;
875     tid = 0;
876
877     if (upper) {
878         apos = inda + *kd;
879         awpos = inda;
880         dpos = apos + *kd;
881         ofdpos = dpos - 1;
882         abdpos = *kd + 1;
883         abofdpos = *kd;
884     } else {
885         apos = inda;
886         awpos = inda + *kd + 1;
887         dpos = apos;
888         ofdpos = dpos + 1;
889         abdpos = 1;
890         abofdpos = 2;
891     }
892
893 /*     Case KD=0: */
894 /*     The matrix is diagonal. We just copy it (convert to "real" for */
895 /*     complex because D is double and the imaginary part should be 0) */
896 /*     and store it in D. A sequential code here is better or */
897 /*     in a parallel environment it might need two cores for D and E */
898
899     if (*kd == 0) {
900         i__1 = *n;
901         for (i__ = 1; i__ <= i__1; ++i__) {
902             i__2 = abdpos + i__ * ab_dim1;
903             d__[i__] = ab[i__2].r;
904 /* L30: */
905         }
906         i__1 = *n - 1;
907         for (i__ = 1; i__ <= i__1; ++i__) {
908             e[i__] = 0.f;
909 /* L40: */
910         }
911
912         hous[1].r = 1.f, hous[1].i = 0.f;
913         work[1].r = 1.f, work[1].i = 0.f;
914         return 0;
915     }
916
917 /*     Case KD=1: */
918 /*     The matrix is already Tridiagonal. We have to make diagonal */
919 /*     and offdiagonal elements real, and store them in D and E. */
920 /*     For that, for real precision just copy the diag and offdiag */
921 /*     to D and E while for the COMPLEX case the bulge chasing is */
922 /*     performed to convert the hermetian tridiagonal to symmetric */
923 /*     tridiagonal. A simpler coversion formula might be used, but then */
924 /*     updating the Q matrix will be required and based if Q is generated */
925 /*     or not this might complicate the story. */
926
927     if (*kd == 1) {
928         i__1 = *n;
929         for (i__ = 1; i__ <= i__1; ++i__) {
930             i__2 = abdpos + i__ * ab_dim1;
931             d__[i__] = ab[i__2].r;
932 /* L50: */
933         }
934
935 /*         make off-diagonal elements real and copy them to E */
936
937         if (upper) {
938             i__1 = *n - 1;
939             for (i__ = 1; i__ <= i__1; ++i__) {
940                 i__2 = abofdpos + (i__ + 1) * ab_dim1;
941                 tmp.r = ab[i__2].r, tmp.i = ab[i__2].i;
942                 abstmp = c_abs(&tmp);
943                 i__2 = abofdpos + (i__ + 1) * ab_dim1;
944                 ab[i__2].r = abstmp, ab[i__2].i = 0.f;
945                 e[i__] = abstmp;
946                 if (abstmp != 0.f) {
947                     q__1.r = tmp.r / abstmp, q__1.i = tmp.i / abstmp;
948                     tmp.r = q__1.r, tmp.i = q__1.i;
949                 } else {
950                     tmp.r = 1.f, tmp.i = 0.f;
951                 }
952                 if (i__ < *n - 1) {
953                     i__2 = abofdpos + (i__ + 2) * ab_dim1;
954                     i__3 = abofdpos + (i__ + 2) * ab_dim1;
955                     q__1.r = ab[i__3].r * tmp.r - ab[i__3].i * tmp.i, q__1.i =
956                              ab[i__3].r * tmp.i + ab[i__3].i * tmp.r;
957                     ab[i__2].r = q__1.r, ab[i__2].i = q__1.i;
958                 }
959 /*                  IF( WANTZ ) THEN */
960 /*                     CALL CSCAL( N, CONJG( TMP ), Q( 1, I+1 ), 1 ) */
961 /*                  END IF */
962 /* L60: */
963             }
964         } else {
965             i__1 = *n - 1;
966             for (i__ = 1; i__ <= i__1; ++i__) {
967                 i__2 = abofdpos + i__ * ab_dim1;
968                 tmp.r = ab[i__2].r, tmp.i = ab[i__2].i;
969                 abstmp = c_abs(&tmp);
970                 i__2 = abofdpos + i__ * ab_dim1;
971                 ab[i__2].r = abstmp, ab[i__2].i = 0.f;
972                 e[i__] = abstmp;
973                 if (abstmp != 0.f) {
974                     q__1.r = tmp.r / abstmp, q__1.i = tmp.i / abstmp;
975                     tmp.r = q__1.r, tmp.i = q__1.i;
976                 } else {
977                     tmp.r = 1.f, tmp.i = 0.f;
978                 }
979                 if (i__ < *n - 1) {
980                     i__2 = abofdpos + (i__ + 1) * ab_dim1;
981                     i__3 = abofdpos + (i__ + 1) * ab_dim1;
982                     q__1.r = ab[i__3].r * tmp.r - ab[i__3].i * tmp.i, q__1.i =
983                              ab[i__3].r * tmp.i + ab[i__3].i * tmp.r;
984                     ab[i__2].r = q__1.r, ab[i__2].i = q__1.i;
985                 }
986 /*                 IF( WANTQ ) THEN */
987 /*                    CALL CSCAL( N, TMP, Q( 1, I+1 ), 1 ) */
988 /*                 END IF */
989 /* L70: */
990             }
991         }
992
993         hous[1].r = 1.f, hous[1].i = 0.f;
994         work[1].r = 1.f, work[1].i = 0.f;
995         return 0;
996     }
997
998 /*     Main code start here. */
999 /*     Reduce the hermitian band of A to a tridiagonal matrix. */
1000
1001     thgrsiz = *n;
1002     grsiz = 1;
1003     shift = 3;
1004     r__1 = (real) (*n) / (real) (*kd) + .5f;
1005     nbtiles = r_int(&r__1);
1006     r__1 = (real) shift / (real) grsiz + .5f;
1007     stepercol = r_int(&r__1);
1008     r__1 = (real) (*n - 1) / (real) thgrsiz + .5f;
1009     thgrnb = r_int(&r__1);
1010
1011     i__1 = *kd + 1;
1012     clacpy_("A", &i__1, n, &ab[ab_offset], ldab, &work[apos], &lda)
1013             ;
1014     claset_("A", kd, n, &c_b1, &c_b1, &work[awpos], &lda);
1015
1016
1017 /*     openMP parallelisation start here */
1018
1019
1020 /*     main bulge chasing loop */
1021
1022     i__1 = thgrnb;
1023     for (thgrid = 1; thgrid <= i__1; ++thgrid) {
1024         stt = (thgrid - 1) * thgrsiz + 1;
1025 /* Computing MIN */
1026         i__2 = stt + thgrsiz - 1, i__3 = *n - 1;
1027         thed = f2cmin(i__2,i__3);
1028         i__2 = *n - 1;
1029         for (i__ = stt; i__ <= i__2; ++i__) {
1030             ed = f2cmin(i__,thed);
1031             if (stt > ed) {
1032                 myexit_();
1033             }
1034             i__3 = stepercol;
1035             for (m = 1; m <= i__3; ++m) {
1036                 st = stt;
1037                 i__4 = ed;
1038                 for (sweepid = st; sweepid <= i__4; ++sweepid) {
1039                     i__5 = grsiz;
1040                     for (k = 1; k <= i__5; ++k) {
1041                         myid = (i__ - sweepid) * (stepercol * grsiz) + (m - 1)
1042                                  * grsiz + k;
1043                         if (myid == 1) {
1044                             ttype = 1;
1045                         } else {
1046                             ttype = myid % 2 + 2;
1047                         }
1048                         if (ttype == 2) {
1049                             colpt = myid / 2 * *kd + sweepid;
1050                             stind = colpt - *kd + 1;
1051                             edind = f2cmin(colpt,*n);
1052                             blklastind = colpt;
1053                         } else {
1054                             colpt = (myid + 1) / 2 * *kd + sweepid;
1055                             stind = colpt - *kd + 1;
1056                             edind = f2cmin(colpt,*n);
1057                             if (stind >= edind - 1 && edind == *n) {
1058                                 blklastind = *n;
1059                             } else {
1060                                 blklastind = 0;
1061                             }
1062                         }
1063
1064 /*                         Call the kernel */
1065
1066                         chb2st_kernels_(uplo, &wantq, &ttype, &stind, &edind,
1067                                  &sweepid, n, kd, &ib, &work[inda], &lda, &
1068                                 hous[indv], &hous[indtau], &ldv, &work[indw + 
1069                                 tid * *kd]);
1070                         if (blklastind >= *n - 1) {
1071                             ++stt;
1072                             myexit_();
1073                         }
1074 /* L140: */
1075                     }
1076 /* L130: */
1077                 }
1078 /* L120: */
1079             }
1080 /* L110: */
1081         }
1082 /* L100: */
1083     }
1084
1085
1086 /*     Copy the diagonal from A to D. Note that D is REAL thus only */
1087 /*     the Real part is needed, the imaginary part should be zero. */
1088
1089     i__1 = *n;
1090     for (i__ = 1; i__ <= i__1; ++i__) {
1091         i__2 = dpos + (i__ - 1) * lda;
1092         d__[i__] = work[i__2].r;
1093 /* L150: */
1094     }
1095
1096 /*     Copy the off diagonal from A to E. Note that E is REAL thus only */
1097 /*     the Real part is needed, the imaginary part should be zero. */
1098
1099     if (upper) {
1100         i__1 = *n - 1;
1101         for (i__ = 1; i__ <= i__1; ++i__) {
1102             i__2 = ofdpos + i__ * lda;
1103             e[i__] = work[i__2].r;
1104 /* L160: */
1105         }
1106     } else {
1107         i__1 = *n - 1;
1108         for (i__ = 1; i__ <= i__1; ++i__) {
1109             i__2 = ofdpos + (i__ - 1) * lda;
1110             e[i__] = work[i__2].r;
1111 /* L170: */
1112         }
1113     }
1114
1115     hous[1].r = (real) lhmin, hous[1].i = 0.f;
1116     work[1].r = (real) lwmin, work[1].i = 0.f;
1117     return 0;
1118
1119 /*     End of CHETRD_HB2ST */
1120
1121 } /* chetrd_hb2st__ */
1122