C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / zheevr_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]/Cd(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__10 = 10;
516 static integer c__1 = 1;
517 static integer c__2 = 2;
518 static integer c__3 = 3;
519 static integer c__4 = 4;
520 static integer c_n1 = -1;
521
522 /* > \brief <b> ZHEEVR_2STAGE computes the eigenvalues and, optionally, the left and/or right eigenvectors for
523  HE matrices</b> */
524
525 /*  @precisions fortran z -> s d c */
526
527 /*  =========== DOCUMENTATION =========== */
528
529 /* Online html documentation available at */
530 /*            http://www.netlib.org/lapack/explore-html/ */
531
532 /* > \htmlonly */
533 /* > Download ZHEEVR_2STAGE + dependencies */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zheevr_
535 2stage.f"> */
536 /* > [TGZ]</a> */
537 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zheevr_
538 2stage.f"> */
539 /* > [ZIP]</a> */
540 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zheevr_
541 2stage.f"> */
542 /* > [TXT]</a> */
543 /* > \endhtmlonly */
544
545 /*  Definition: */
546 /*  =========== */
547
548 /*       SUBROUTINE ZHEEVR_2STAGE( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, */
549 /*                                 IL, IU, ABSTOL, M, W, Z, LDZ, ISUPPZ, */
550 /*                                 WORK, LWORK, RWORK, LRWORK, IWORK, */
551 /*                                 LIWORK, INFO ) */
552
553 /*       IMPLICIT NONE */
554
555 /*       CHARACTER          JOBZ, RANGE, UPLO */
556 /*       INTEGER            IL, INFO, IU, LDA, LDZ, LIWORK, LRWORK, LWORK, */
557 /*      $                   M, N */
558 /*       DOUBLE PRECISION   ABSTOL, VL, VU */
559 /*       INTEGER            ISUPPZ( * ), IWORK( * ) */
560 /*       DOUBLE PRECISION   RWORK( * ), W( * ) */
561 /*       COMPLEX*16         A( LDA, * ), WORK( * ), Z( LDZ, * ) */
562
563
564 /* > \par Purpose: */
565 /*  ============= */
566 /* > */
567 /* > \verbatim */
568 /* > */
569 /* > ZHEEVR_2STAGE computes selected eigenvalues and, optionally, eigenvectors */
570 /* > of a complex Hermitian matrix A using the 2stage technique for */
571 /* > the reduction to tridiagonal.  Eigenvalues and eigenvectors can */
572 /* > be selected by specifying either a range of values or a range of */
573 /* > indices for the desired eigenvalues. */
574 /* > */
575 /* > ZHEEVR_2STAGE first reduces the matrix A to tridiagonal form T with a call */
576 /* > to ZHETRD.  Then, whenever possible, ZHEEVR_2STAGE calls ZSTEMR to compute */
577 /* > eigenspectrum using Relatively Robust Representations.  ZSTEMR */
578 /* > computes eigenvalues by the dqds algorithm, while orthogonal */
579 /* > eigenvectors are computed from various "good" L D L^T representations */
580 /* > (also known as Relatively Robust Representations). Gram-Schmidt */
581 /* > orthogonalization is avoided as far as possible. More specifically, */
582 /* > the various steps of the algorithm are as follows. */
583 /* > */
584 /* > For each unreduced block (submatrix) of T, */
585 /* >    (a) Compute T - sigma I  = L D L^T, so that L and D */
586 /* >        define all the wanted eigenvalues to high relative accuracy. */
587 /* >        This means that small relative changes in the entries of D and L */
588 /* >        cause only small relative changes in the eigenvalues and */
589 /* >        eigenvectors. The standard (unfactored) representation of the */
590 /* >        tridiagonal matrix T does not have this property in general. */
591 /* >    (b) Compute the eigenvalues to suitable accuracy. */
592 /* >        If the eigenvectors are desired, the algorithm attains full */
593 /* >        accuracy of the computed eigenvalues only right before */
594 /* >        the corresponding vectors have to be computed, see steps c) and d). */
595 /* >    (c) For each cluster of close eigenvalues, select a new */
596 /* >        shift close to the cluster, find a new factorization, and refine */
597 /* >        the shifted eigenvalues to suitable accuracy. */
598 /* >    (d) For each eigenvalue with a large enough relative separation compute */
599 /* >        the corresponding eigenvector by forming a rank revealing twisted */
600 /* >        factorization. Go back to (c) for any clusters that remain. */
601 /* > */
602 /* > The desired accuracy of the output can be specified by the input */
603 /* > parameter ABSTOL. */
604 /* > */
605 /* > For more details, see DSTEMR's documentation and: */
606 /* > - Inderjit S. Dhillon and Beresford N. Parlett: "Multiple representations */
607 /* >   to compute orthogonal eigenvectors of symmetric tridiagonal matrices," */
608 /* >   Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004. */
609 /* > - Inderjit Dhillon and Beresford Parlett: "Orthogonal Eigenvectors and */
610 /* >   Relative Gaps," SIAM Journal on Matrix Analysis and Applications, Vol. 25, */
611 /* >   2004.  Also LAPACK Working Note 154. */
612 /* > - Inderjit Dhillon: "A new O(n^2) algorithm for the symmetric */
613 /* >   tridiagonal eigenvalue/eigenvector problem", */
614 /* >   Computer Science Division Technical Report No. UCB/CSD-97-971, */
615 /* >   UC Berkeley, May 1997. */
616 /* > */
617 /* > */
618 /* > Note 1 : ZHEEVR_2STAGE calls ZSTEMR when the full spectrum is requested */
619 /* > on machines which conform to the ieee-754 floating point standard. */
620 /* > ZHEEVR_2STAGE calls DSTEBZ and ZSTEIN on non-ieee machines and */
621 /* > when partial spectrum requests are made. */
622 /* > */
623 /* > Normal execution of ZSTEMR may create NaNs and infinities and */
624 /* > hence may abort due to a floating point exception in environments */
625 /* > which do not handle NaNs and infinities in the ieee standard default */
626 /* > manner. */
627 /* > \endverbatim */
628
629 /*  Arguments: */
630 /*  ========== */
631
632 /* > \param[in] JOBZ */
633 /* > \verbatim */
634 /* >          JOBZ is CHARACTER*1 */
635 /* >          = 'N':  Compute eigenvalues only; */
636 /* >          = 'V':  Compute eigenvalues and eigenvectors. */
637 /* >                  Not available in this release. */
638 /* > \endverbatim */
639 /* > */
640 /* > \param[in] RANGE */
641 /* > \verbatim */
642 /* >          RANGE is CHARACTER*1 */
643 /* >          = 'A': all eigenvalues will be found. */
644 /* >          = 'V': all eigenvalues in the half-open interval (VL,VU] */
645 /* >                 will be found. */
646 /* >          = 'I': the IL-th through IU-th eigenvalues will be found. */
647 /* >          For RANGE = 'V' or 'I' and IU - IL < N - 1, DSTEBZ and */
648 /* >          ZSTEIN are called */
649 /* > \endverbatim */
650 /* > */
651 /* > \param[in] UPLO */
652 /* > \verbatim */
653 /* >          UPLO is CHARACTER*1 */
654 /* >          = 'U':  Upper triangle of A is stored; */
655 /* >          = 'L':  Lower triangle of A is stored. */
656 /* > \endverbatim */
657 /* > */
658 /* > \param[in] N */
659 /* > \verbatim */
660 /* >          N is INTEGER */
661 /* >          The order of the matrix A.  N >= 0. */
662 /* > \endverbatim */
663 /* > */
664 /* > \param[in,out] A */
665 /* > \verbatim */
666 /* >          A is COMPLEX*16 array, dimension (LDA, N) */
667 /* >          On entry, the Hermitian matrix A.  If UPLO = 'U', the */
668 /* >          leading N-by-N upper triangular part of A contains the */
669 /* >          upper triangular part of the matrix A.  If UPLO = 'L', */
670 /* >          the leading N-by-N lower triangular part of A contains */
671 /* >          the lower triangular part of the matrix A. */
672 /* >          On exit, the lower triangle (if UPLO='L') or the upper */
673 /* >          triangle (if UPLO='U') of A, including the diagonal, is */
674 /* >          destroyed. */
675 /* > \endverbatim */
676 /* > */
677 /* > \param[in] LDA */
678 /* > \verbatim */
679 /* >          LDA is INTEGER */
680 /* >          The leading dimension of the array A.  LDA >= f2cmax(1,N). */
681 /* > \endverbatim */
682 /* > */
683 /* > \param[in] VL */
684 /* > \verbatim */
685 /* >          VL is DOUBLE PRECISION */
686 /* >          If RANGE='V', the lower bound of the interval to */
687 /* >          be searched for eigenvalues. VL < VU. */
688 /* >          Not referenced if RANGE = 'A' or 'I'. */
689 /* > \endverbatim */
690 /* > */
691 /* > \param[in] VU */
692 /* > \verbatim */
693 /* >          VU is DOUBLE PRECISION */
694 /* >          If RANGE='V', the upper bound of the interval to */
695 /* >          be searched for eigenvalues. VL < VU. */
696 /* >          Not referenced if RANGE = 'A' or 'I'. */
697 /* > \endverbatim */
698 /* > */
699 /* > \param[in] IL */
700 /* > \verbatim */
701 /* >          IL is INTEGER */
702 /* >          If RANGE='I', the index of the */
703 /* >          smallest eigenvalue to be returned. */
704 /* >          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. */
705 /* >          Not referenced if RANGE = 'A' or 'V'. */
706 /* > \endverbatim */
707 /* > */
708 /* > \param[in] IU */
709 /* > \verbatim */
710 /* >          IU is INTEGER */
711 /* >          If RANGE='I', the index of the */
712 /* >          largest eigenvalue to be returned. */
713 /* >          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. */
714 /* >          Not referenced if RANGE = 'A' or 'V'. */
715 /* > \endverbatim */
716 /* > */
717 /* > \param[in] ABSTOL */
718 /* > \verbatim */
719 /* >          ABSTOL is DOUBLE PRECISION */
720 /* >          The absolute error tolerance for the eigenvalues. */
721 /* >          An approximate eigenvalue is accepted as converged */
722 /* >          when it is determined to lie in an interval [a,b] */
723 /* >          of width less than or equal to */
724 /* > */
725 /* >                  ABSTOL + EPS *   f2cmax( |a|,|b| ) , */
726 /* > */
727 /* >          where EPS is the machine precision.  If ABSTOL is less than */
728 /* >          or equal to zero, then  EPS*|T|  will be used in its place, */
729 /* >          where |T| is the 1-norm of the tridiagonal matrix obtained */
730 /* >          by reducing A to tridiagonal form. */
731 /* > */
732 /* >          See "Computing Small Singular Values of Bidiagonal Matrices */
733 /* >          with Guaranteed High Relative Accuracy," by Demmel and */
734 /* >          Kahan, LAPACK Working Note #3. */
735 /* > */
736 /* >          If high relative accuracy is important, set ABSTOL to */
737 /* >          DLAMCH( 'Safe minimum' ).  Doing so will guarantee that */
738 /* >          eigenvalues are computed to high relative accuracy when */
739 /* >          possible in future releases.  The current code does not */
740 /* >          make any guarantees about high relative accuracy, but */
741 /* >          future releases will. See J. Barlow and J. Demmel, */
742 /* >          "Computing Accurate Eigensystems of Scaled Diagonally */
743 /* >          Dominant Matrices", LAPACK Working Note #7, for a discussion */
744 /* >          of which matrices define their eigenvalues to high relative */
745 /* >          accuracy. */
746 /* > \endverbatim */
747 /* > */
748 /* > \param[out] M */
749 /* > \verbatim */
750 /* >          M is INTEGER */
751 /* >          The total number of eigenvalues found.  0 <= M <= N. */
752 /* >          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1. */
753 /* > \endverbatim */
754 /* > */
755 /* > \param[out] W */
756 /* > \verbatim */
757 /* >          W is DOUBLE PRECISION array, dimension (N) */
758 /* >          The first M elements contain the selected eigenvalues in */
759 /* >          ascending order. */
760 /* > \endverbatim */
761 /* > */
762 /* > \param[out] Z */
763 /* > \verbatim */
764 /* >          Z is COMPLEX*16 array, dimension (LDZ, f2cmax(1,M)) */
765 /* >          If JOBZ = 'V', then if INFO = 0, the first M columns of Z */
766 /* >          contain the orthonormal eigenvectors of the matrix A */
767 /* >          corresponding to the selected eigenvalues, with the i-th */
768 /* >          column of Z holding the eigenvector associated with W(i). */
769 /* >          If JOBZ = 'N', then Z is not referenced. */
770 /* >          Note: the user must ensure that at least f2cmax(1,M) columns are */
771 /* >          supplied in the array Z; if RANGE = 'V', the exact value of M */
772 /* >          is not known in advance and an upper bound must be used. */
773 /* > \endverbatim */
774 /* > */
775 /* > \param[in] LDZ */
776 /* > \verbatim */
777 /* >          LDZ is INTEGER */
778 /* >          The leading dimension of the array Z.  LDZ >= 1, and if */
779 /* >          JOBZ = 'V', LDZ >= f2cmax(1,N). */
780 /* > \endverbatim */
781 /* > */
782 /* > \param[out] ISUPPZ */
783 /* > \verbatim */
784 /* >          ISUPPZ is INTEGER array, dimension ( 2*f2cmax(1,M) ) */
785 /* >          The support of the eigenvectors in Z, i.e., the indices */
786 /* >          indicating the nonzero elements in Z. The i-th eigenvector */
787 /* >          is nonzero only in elements ISUPPZ( 2*i-1 ) through */
788 /* >          ISUPPZ( 2*i ). This is an output of ZSTEMR (tridiagonal */
789 /* >          matrix). The support of the eigenvectors of A is typically */
790 /* >          1:N because of the unitary transformations applied by ZUNMTR. */
791 /* >          Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1 */
792 /* > \endverbatim */
793 /* > */
794 /* > \param[out] WORK */
795 /* > \verbatim */
796 /* >          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) */
797 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
798 /* > \endverbatim */
799 /* > */
800 /* > \param[in] LWORK */
801 /* > \verbatim */
802 /* >          LWORK is INTEGER */
803 /* >          The dimension of the array WORK. */
804 /* >          If JOBZ = 'N' and N > 1, LWORK must be queried. */
805 /* >                                   LWORK = MAX(1, 26*N, dimension) where */
806 /* >                                   dimension = f2cmax(stage1,stage2) + (KD+1)*N + N */
807 /* >                                             = N*KD + N*f2cmax(KD+1,FACTOPTNB) */
808 /* >                                               + f2cmax(2*KD*KD, KD*NTHREADS) */
809 /* >                                               + (KD+1)*N + N */
810 /* >                                   where KD is the blocking size of the reduction, */
811 /* >                                   FACTOPTNB is the blocking used by the QR or LQ */
812 /* >                                   algorithm, usually FACTOPTNB=128 is a good choice */
813 /* >                                   NTHREADS is the number of threads used when */
814 /* >                                   openMP compilation is enabled, otherwise =1. */
815 /* >          If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available */
816 /* > */
817 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
818 /* >          only calculates the optimal sizes of the WORK, RWORK and */
819 /* >          IWORK arrays, returns these values as the first entries of */
820 /* >          the WORK, RWORK and IWORK arrays, and no error message */
821 /* >          related to LWORK or LRWORK or LIWORK is issued by XERBLA. */
822 /* > \endverbatim */
823 /* > */
824 /* > \param[out] RWORK */
825 /* > \verbatim */
826 /* >          RWORK is DOUBLE PRECISION array, dimension (MAX(1,LRWORK)) */
827 /* >          On exit, if INFO = 0, RWORK(1) returns the optimal */
828 /* >          (and minimal) LRWORK. */
829 /* > \endverbatim */
830 /* > */
831 /* > \param[in] LRWORK */
832 /* > \verbatim */
833 /* >          LRWORK is INTEGER */
834 /* >          The length of the array RWORK.  LRWORK >= f2cmax(1,24*N). */
835 /* > */
836 /* >          If LRWORK = -1, then a workspace query is assumed; the */
837 /* >          routine only calculates the optimal sizes of the WORK, RWORK */
838 /* >          and IWORK arrays, returns these values as the first entries */
839 /* >          of the WORK, RWORK and IWORK arrays, and no error message */
840 /* >          related to LWORK or LRWORK or LIWORK is issued by XERBLA. */
841 /* > \endverbatim */
842 /* > */
843 /* > \param[out] IWORK */
844 /* > \verbatim */
845 /* >          IWORK is INTEGER array, dimension (MAX(1,LIWORK)) */
846 /* >          On exit, if INFO = 0, IWORK(1) returns the optimal */
847 /* >          (and minimal) LIWORK. */
848 /* > \endverbatim */
849 /* > */
850 /* > \param[in] LIWORK */
851 /* > \verbatim */
852 /* >          LIWORK is INTEGER */
853 /* >          The dimension of the array IWORK.  LIWORK >= f2cmax(1,10*N). */
854 /* > */
855 /* >          If LIWORK = -1, then a workspace query is assumed; the */
856 /* >          routine only calculates the optimal sizes of the WORK, RWORK */
857 /* >          and IWORK arrays, returns these values as the first entries */
858 /* >          of the WORK, RWORK and IWORK arrays, and no error message */
859 /* >          related to LWORK or LRWORK or LIWORK is issued by XERBLA. */
860 /* > \endverbatim */
861 /* > */
862 /* > \param[out] INFO */
863 /* > \verbatim */
864 /* >          INFO is INTEGER */
865 /* >          = 0:  successful exit */
866 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value */
867 /* >          > 0:  Internal error */
868 /* > \endverbatim */
869
870 /*  Authors: */
871 /*  ======== */
872
873 /* > \author Univ. of Tennessee */
874 /* > \author Univ. of California Berkeley */
875 /* > \author Univ. of Colorado Denver */
876 /* > \author NAG Ltd. */
877
878 /* > \date June 2016 */
879
880 /* > \ingroup complex16HEeigen */
881
882 /* > \par Contributors: */
883 /*  ================== */
884 /* > */
885 /* >     Inderjit Dhillon, IBM Almaden, USA \n */
886 /* >     Osni Marques, LBNL/NERSC, USA \n */
887 /* >     Ken Stanley, Computer Science Division, University of */
888 /* >       California at Berkeley, USA \n */
889 /* >     Jason Riedy, Computer Science Division, University of */
890 /* >       California at Berkeley, USA \n */
891 /* > */
892 /* > \par Further Details: */
893 /*  ===================== */
894 /* > */
895 /* > \verbatim */
896 /* > */
897 /* >  All details about the 2stage techniques are available in: */
898 /* > */
899 /* >  Azzam Haidar, Hatem Ltaief, and Jack Dongarra. */
900 /* >  Parallel reduction to condensed forms for symmetric eigenvalue problems */
901 /* >  using aggregated fine-grained and memory-aware kernels. In Proceedings */
902 /* >  of 2011 International Conference for High Performance Computing, */
903 /* >  Networking, Storage and Analysis (SC '11), New York, NY, USA, */
904 /* >  Article 8 , 11 pages. */
905 /* >  http://doi.acm.org/10.1145/2063384.2063394 */
906 /* > */
907 /* >  A. Haidar, J. Kurzak, P. Luszczek, 2013. */
908 /* >  An improved parallel singular value algorithm and its implementation */
909 /* >  for multicore hardware, In Proceedings of 2013 International Conference */
910 /* >  for High Performance Computing, Networking, Storage and Analysis (SC '13). */
911 /* >  Denver, Colorado, USA, 2013. */
912 /* >  Article 90, 12 pages. */
913 /* >  http://doi.acm.org/10.1145/2503210.2503292 */
914 /* > */
915 /* >  A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra. */
916 /* >  A novel hybrid CPU-GPU generalized eigensolver for electronic structure */
917 /* >  calculations based on fine-grained memory aware tasks. */
918 /* >  International Journal of High Performance Computing Applications. */
919 /* >  Volume 28 Issue 2, Pages 196-209, May 2014. */
920 /* >  http://hpc.sagepub.com/content/28/2/196 */
921 /* > */
922 /* > \endverbatim */
923
924 /*  ===================================================================== */
925 /* Subroutine */ int zheevr_2stage_(char *jobz, char *range, char *uplo, 
926         integer *n, doublecomplex *a, integer *lda, doublereal *vl, 
927         doublereal *vu, integer *il, integer *iu, doublereal *abstol, integer 
928         *m, doublereal *w, doublecomplex *z__, integer *ldz, integer *isuppz, 
929         doublecomplex *work, integer *lwork, doublereal *rwork, integer *
930         lrwork, integer *iwork, integer *liwork, integer *info)
931 {
932     /* System generated locals */
933     integer a_dim1, a_offset, z_dim1, z_offset, i__1, i__2;
934     doublereal d__1, d__2;
935
936     /* Local variables */
937     extern integer ilaenv2stage_(integer *, char *, char *, integer *, 
938             integer *, integer *, integer *);
939     doublereal anrm;
940     integer imax;
941     doublereal rmin, rmax;
942     logical test;
943     integer itmp1, i__, j;
944     extern /* Subroutine */ int dscal_(integer *, doublereal *, doublereal *, 
945             integer *);
946     integer indrd, indre;
947     doublereal sigma;
948     extern logical lsame_(char *, char *);
949     integer iinfo;
950     extern /* Subroutine */ int zhetrd_2stage_(char *, char *, integer *, 
951             doublecomplex *, integer *, doublereal *, doublereal *, 
952             doublecomplex *, doublecomplex *, integer *, doublecomplex *, 
953             integer *, integer *);
954     char order[1];
955     integer indwk, lhtrd;
956     extern /* Subroutine */ int dcopy_(integer *, doublereal *, integer *, 
957             doublereal *, integer *);
958     integer lwmin;
959     logical lower;
960     integer lwtrd;
961     logical wantz;
962     extern /* Subroutine */ int zswap_(integer *, doublecomplex *, integer *, 
963             doublecomplex *, integer *);
964     integer ib, kd, jj;
965     extern doublereal dlamch_(char *);
966     logical alleig, indeig;
967     integer iscale, ieeeok, indibl, indrdd, indifl, indree;
968     logical valeig;
969     doublereal safmin;
970     extern integer ilaenv_(integer *, char *, char *, integer *, integer *, 
971             integer *, integer *, ftnlen, ftnlen);
972     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen), zdscal_(
973             integer *, doublereal *, doublecomplex *, integer *);
974     doublereal abstll, bignum;
975     integer indtau, indisp;
976     extern /* Subroutine */ int dsterf_(integer *, doublereal *, doublereal *,
977              integer *);
978     integer indiwo, indwkn;
979     extern /* Subroutine */ int dstebz_(char *, char *, integer *, doublereal 
980             *, doublereal *, integer *, integer *, doublereal *, doublereal *,
981              doublereal *, integer *, integer *, doublereal *, integer *, 
982             integer *, doublereal *, integer *, integer *);
983     integer indrwk, liwmin;
984     logical tryrac;
985     integer lrwmin, llwrkn, llwork, nsplit;
986     doublereal smlnum;
987     extern /* Subroutine */ int zstein_(integer *, doublereal *, doublereal *,
988              integer *, doublereal *, integer *, integer *, doublecomplex *, 
989             integer *, doublereal *, integer *, integer *, integer *);
990     logical lquery;
991     extern doublereal zlansy_(char *, char *, integer *, doublecomplex *, 
992             integer *, doublereal *);
993     extern /* Subroutine */ int zstemr_(char *, char *, integer *, doublereal 
994             *, doublereal *, doublereal *, doublereal *, integer *, integer *,
995              integer *, doublereal *, doublecomplex *, integer *, integer *, 
996             integer *, logical *, doublereal *, integer *, integer *, integer 
997             *, integer *), zunmtr_(char *, char *, char *, 
998             integer *, integer *, doublecomplex *, integer *, doublecomplex *,
999              doublecomplex *, integer *, doublecomplex *, integer *, integer *
1000             );
1001     doublereal eps, vll, vuu;
1002     integer indhous, llrwork;
1003     doublereal tmp1;
1004
1005
1006
1007 /*  -- LAPACK driver routine (version 3.8.0) -- */
1008 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
1009 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
1010 /*     June 2016 */
1011
1012
1013 /*  ===================================================================== */
1014
1015
1016 /*     Test the input parameters. */
1017
1018     /* Parameter adjustments */
1019     a_dim1 = *lda;
1020     a_offset = 1 + a_dim1 * 1;
1021     a -= a_offset;
1022     --w;
1023     z_dim1 = *ldz;
1024     z_offset = 1 + z_dim1 * 1;
1025     z__ -= z_offset;
1026     --isuppz;
1027     --work;
1028     --rwork;
1029     --iwork;
1030
1031     /* Function Body */
1032     ieeeok = ilaenv_(&c__10, "ZHEEVR", "N", &c__1, &c__2, &c__3, &c__4, (
1033             ftnlen)6, (ftnlen)1);
1034
1035     lower = lsame_(uplo, "L");
1036     wantz = lsame_(jobz, "V");
1037     alleig = lsame_(range, "A");
1038     valeig = lsame_(range, "V");
1039     indeig = lsame_(range, "I");
1040
1041     lquery = *lwork == -1 || *lrwork == -1 || *liwork == -1;
1042
1043     kd = ilaenv2stage_(&c__1, "ZHETRD_2STAGE", jobz, n, &c_n1, &c_n1, &c_n1);
1044     ib = ilaenv2stage_(&c__2, "ZHETRD_2STAGE", jobz, n, &kd, &c_n1, &c_n1);
1045     lhtrd = ilaenv2stage_(&c__3, "ZHETRD_2STAGE", jobz, n, &kd, &ib, &c_n1);
1046     lwtrd = ilaenv2stage_(&c__4, "ZHETRD_2STAGE", jobz, n, &kd, &ib, &c_n1);
1047     lwmin = *n + lhtrd + lwtrd;
1048 /* Computing MAX */
1049     i__1 = 1, i__2 = *n * 24;
1050     lrwmin = f2cmax(i__1,i__2);
1051 /* Computing MAX */
1052     i__1 = 1, i__2 = *n * 10;
1053     liwmin = f2cmax(i__1,i__2);
1054
1055     *info = 0;
1056     if (! lsame_(jobz, "N")) {
1057         *info = -1;
1058     } else if (! (alleig || valeig || indeig)) {
1059         *info = -2;
1060     } else if (! (lower || lsame_(uplo, "U"))) {
1061         *info = -3;
1062     } else if (*n < 0) {
1063         *info = -4;
1064     } else if (*lda < f2cmax(1,*n)) {
1065         *info = -6;
1066     } else {
1067         if (valeig) {
1068             if (*n > 0 && *vu <= *vl) {
1069                 *info = -8;
1070             }
1071         } else if (indeig) {
1072             if (*il < 1 || *il > f2cmax(1,*n)) {
1073                 *info = -9;
1074             } else if (*iu < f2cmin(*n,*il) || *iu > *n) {
1075                 *info = -10;
1076             }
1077         }
1078     }
1079     if (*info == 0) {
1080         if (*ldz < 1 || wantz && *ldz < *n) {
1081             *info = -15;
1082         }
1083     }
1084
1085     if (*info == 0) {
1086         work[1].r = (doublereal) lwmin, work[1].i = 0.;
1087         rwork[1] = (doublereal) lrwmin;
1088         iwork[1] = liwmin;
1089
1090         if (*lwork < lwmin && ! lquery) {
1091             *info = -18;
1092         } else if (*lrwork < lrwmin && ! lquery) {
1093             *info = -20;
1094         } else if (*liwork < liwmin && ! lquery) {
1095             *info = -22;
1096         }
1097     }
1098
1099     if (*info != 0) {
1100         i__1 = -(*info);
1101         xerbla_("ZHEEVR_2STAGE", &i__1, (ftnlen)13);
1102         return 0;
1103     } else if (lquery) {
1104         return 0;
1105     }
1106
1107 /*     Quick return if possible */
1108
1109     *m = 0;
1110     if (*n == 0) {
1111         work[1].r = 1., work[1].i = 0.;
1112         return 0;
1113     }
1114
1115     if (*n == 1) {
1116         work[1].r = 2., work[1].i = 0.;
1117         if (alleig || indeig) {
1118             *m = 1;
1119             i__1 = a_dim1 + 1;
1120             w[1] = a[i__1].r;
1121         } else {
1122             i__1 = a_dim1 + 1;
1123             i__2 = a_dim1 + 1;
1124             if (*vl < a[i__1].r && *vu >= a[i__2].r) {
1125                 *m = 1;
1126                 i__1 = a_dim1 + 1;
1127                 w[1] = a[i__1].r;
1128             }
1129         }
1130         if (wantz) {
1131             i__1 = z_dim1 + 1;
1132             z__[i__1].r = 1., z__[i__1].i = 0.;
1133             isuppz[1] = 1;
1134             isuppz[2] = 1;
1135         }
1136         return 0;
1137     }
1138
1139 /*     Get machine constants. */
1140
1141     safmin = dlamch_("Safe minimum");
1142     eps = dlamch_("Precision");
1143     smlnum = safmin / eps;
1144     bignum = 1. / smlnum;
1145     rmin = sqrt(smlnum);
1146 /* Computing MIN */
1147     d__1 = sqrt(bignum), d__2 = 1. / sqrt(sqrt(safmin));
1148     rmax = f2cmin(d__1,d__2);
1149
1150 /*     Scale matrix to allowable range, if necessary. */
1151
1152     iscale = 0;
1153     abstll = *abstol;
1154     if (valeig) {
1155         vll = *vl;
1156         vuu = *vu;
1157     }
1158     anrm = zlansy_("M", uplo, n, &a[a_offset], lda, &rwork[1]);
1159     if (anrm > 0. && anrm < rmin) {
1160         iscale = 1;
1161         sigma = rmin / anrm;
1162     } else if (anrm > rmax) {
1163         iscale = 1;
1164         sigma = rmax / anrm;
1165     }
1166     if (iscale == 1) {
1167         if (lower) {
1168             i__1 = *n;
1169             for (j = 1; j <= i__1; ++j) {
1170                 i__2 = *n - j + 1;
1171                 zdscal_(&i__2, &sigma, &a[j + j * a_dim1], &c__1);
1172 /* L10: */
1173             }
1174         } else {
1175             i__1 = *n;
1176             for (j = 1; j <= i__1; ++j) {
1177                 zdscal_(&j, &sigma, &a[j * a_dim1 + 1], &c__1);
1178 /* L20: */
1179             }
1180         }
1181         if (*abstol > 0.) {
1182             abstll = *abstol * sigma;
1183         }
1184         if (valeig) {
1185             vll = *vl * sigma;
1186             vuu = *vu * sigma;
1187         }
1188     }
1189 /*     Initialize indices into workspaces.  Note: The IWORK indices are */
1190 /*     used only if DSTERF or ZSTEMR fail. */
1191 /*     WORK(INDTAU:INDTAU+N-1) stores the complex scalar factors of the */
1192 /*     elementary reflectors used in ZHETRD. */
1193     indtau = 1;
1194 /*     INDWK is the starting offset of the remaining complex workspace, */
1195 /*     and LLWORK is the remaining complex workspace size. */
1196     indhous = indtau + *n;
1197     indwk = indhous + lhtrd;
1198     llwork = *lwork - indwk + 1;
1199 /*     RWORK(INDRD:INDRD+N-1) stores the real tridiagonal's diagonal */
1200 /*     entries. */
1201     indrd = 1;
1202 /*     RWORK(INDRE:INDRE+N-1) stores the off-diagonal entries of the */
1203 /*     tridiagonal matrix from ZHETRD. */
1204     indre = indrd + *n;
1205 /*     RWORK(INDRDD:INDRDD+N-1) is a copy of the diagonal entries over */
1206 /*     -written by ZSTEMR (the DSTERF path copies the diagonal to W). */
1207     indrdd = indre + *n;
1208 /*     RWORK(INDREE:INDREE+N-1) is a copy of the off-diagonal entries over */
1209 /*     -written while computing the eigenvalues in DSTERF and ZSTEMR. */
1210     indree = indrdd + *n;
1211 /*     INDRWK is the starting offset of the left-over real workspace, and */
1212 /*     LLRWORK is the remaining workspace size. */
1213     indrwk = indree + *n;
1214     llrwork = *lrwork - indrwk + 1;
1215 /*     IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in DSTEBZ and */
1216 /*     stores the block indices of each of the M<=N eigenvalues. */
1217     indibl = 1;
1218 /*     IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in DSTEBZ and */
1219 /*     stores the starting and finishing indices of each block. */
1220     indisp = indibl + *n;
1221 /*     IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors */
1222 /*     that corresponding to eigenvectors that fail to converge in */
1223 /*     ZSTEIN.  This information is discarded; if any fail, the driver */
1224 /*     returns INFO > 0. */
1225     indifl = indisp + *n;
1226 /*     INDIWO is the offset of the remaining integer workspace. */
1227     indiwo = indifl + *n;
1228
1229 /*     Call ZHETRD_2STAGE to reduce Hermitian matrix to tridiagonal form. */
1230
1231     zhetrd_2stage_(jobz, uplo, n, &a[a_offset], lda, &rwork[indrd], &rwork[
1232             indre], &work[indtau], &work[indhous], &lhtrd, &work[indwk], &
1233             llwork, &iinfo);
1234
1235 /*     If all eigenvalues are desired */
1236 /*     then call DSTERF or ZSTEMR and ZUNMTR. */
1237
1238     test = FALSE_;
1239     if (indeig) {
1240         if (*il == 1 && *iu == *n) {
1241             test = TRUE_;
1242         }
1243     }
1244     if ((alleig || test) && ieeeok == 1) {
1245         if (! wantz) {
1246             dcopy_(n, &rwork[indrd], &c__1, &w[1], &c__1);
1247             i__1 = *n - 1;
1248             dcopy_(&i__1, &rwork[indre], &c__1, &rwork[indree], &c__1);
1249             dsterf_(n, &w[1], &rwork[indree], info);
1250         } else {
1251             i__1 = *n - 1;
1252             dcopy_(&i__1, &rwork[indre], &c__1, &rwork[indree], &c__1);
1253             dcopy_(n, &rwork[indrd], &c__1, &rwork[indrdd], &c__1);
1254
1255             if (*abstol <= *n * 2. * eps) {
1256                 tryrac = TRUE_;
1257             } else {
1258                 tryrac = FALSE_;
1259             }
1260             zstemr_(jobz, "A", n, &rwork[indrdd], &rwork[indree], vl, vu, il, 
1261                     iu, m, &w[1], &z__[z_offset], ldz, n, &isuppz[1], &tryrac,
1262                      &rwork[indrwk], &llrwork, &iwork[1], liwork, info);
1263
1264 /*           Apply unitary matrix used in reduction to tridiagonal */
1265 /*           form to eigenvectors returned by ZSTEMR. */
1266
1267             if (wantz && *info == 0) {
1268                 indwkn = indwk;
1269                 llwrkn = *lwork - indwkn + 1;
1270                 zunmtr_("L", uplo, "N", n, m, &a[a_offset], lda, &work[indtau]
1271                         , &z__[z_offset], ldz, &work[indwkn], &llwrkn, &iinfo);
1272             }
1273         }
1274
1275
1276         if (*info == 0) {
1277             *m = *n;
1278             goto L30;
1279         }
1280         *info = 0;
1281     }
1282
1283 /*     Otherwise, call DSTEBZ and, if eigenvectors are desired, ZSTEIN. */
1284 /*     Also call DSTEBZ and ZSTEIN if ZSTEMR fails. */
1285
1286     if (wantz) {
1287         *(unsigned char *)order = 'B';
1288     } else {
1289         *(unsigned char *)order = 'E';
1290     }
1291     dstebz_(range, order, n, &vll, &vuu, il, iu, &abstll, &rwork[indrd], &
1292             rwork[indre], m, &nsplit, &w[1], &iwork[indibl], &iwork[indisp], &
1293             rwork[indrwk], &iwork[indiwo], info);
1294
1295     if (wantz) {
1296         zstein_(n, &rwork[indrd], &rwork[indre], m, &w[1], &iwork[indibl], &
1297                 iwork[indisp], &z__[z_offset], ldz, &rwork[indrwk], &iwork[
1298                 indiwo], &iwork[indifl], info);
1299
1300 /*        Apply unitary matrix used in reduction to tridiagonal */
1301 /*        form to eigenvectors returned by ZSTEIN. */
1302
1303         indwkn = indwk;
1304         llwrkn = *lwork - indwkn + 1;
1305         zunmtr_("L", uplo, "N", n, m, &a[a_offset], lda, &work[indtau], &z__[
1306                 z_offset], ldz, &work[indwkn], &llwrkn, &iinfo);
1307     }
1308
1309 /*     If matrix was scaled, then rescale eigenvalues appropriately. */
1310
1311 L30:
1312     if (iscale == 1) {
1313         if (*info == 0) {
1314             imax = *m;
1315         } else {
1316             imax = *info - 1;
1317         }
1318         d__1 = 1. / sigma;
1319         dscal_(&imax, &d__1, &w[1], &c__1);
1320     }
1321
1322 /*     If eigenvalues are not in order, then sort them, along with */
1323 /*     eigenvectors. */
1324
1325     if (wantz) {
1326         i__1 = *m - 1;
1327         for (j = 1; j <= i__1; ++j) {
1328             i__ = 0;
1329             tmp1 = w[j];
1330             i__2 = *m;
1331             for (jj = j + 1; jj <= i__2; ++jj) {
1332                 if (w[jj] < tmp1) {
1333                     i__ = jj;
1334                     tmp1 = w[jj];
1335                 }
1336 /* L40: */
1337             }
1338
1339             if (i__ != 0) {
1340                 itmp1 = iwork[indibl + i__ - 1];
1341                 w[i__] = w[j];
1342                 iwork[indibl + i__ - 1] = iwork[indibl + j - 1];
1343                 w[j] = tmp1;
1344                 iwork[indibl + j - 1] = itmp1;
1345                 zswap_(n, &z__[i__ * z_dim1 + 1], &c__1, &z__[j * z_dim1 + 1],
1346                          &c__1);
1347             }
1348 /* L50: */
1349         }
1350     }
1351
1352 /*     Set WORK(1) to optimal workspace size. */
1353
1354     work[1].r = (doublereal) lwmin, work[1].i = 0.;
1355     rwork[1] = (doublereal) lrwmin;
1356     iwork[1] = liwmin;
1357
1358     return 0;
1359
1360 /*     End of ZHEEVR_2STAGE */
1361
1362 } /* zheevr_2stage__ */
1363