C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / cheevr.c
1 #include <math.h>
2 #include <stdlib.h>
3 #include <string.h>
4 #include <stdio.h>
5 #include <complex.h>
6 #ifdef complex
7 #undef complex
8 #endif
9 #ifdef I
10 #undef I
11 #endif
12
13 #if defined(_WIN64)
14 typedef long long BLASLONG;
15 typedef unsigned long long BLASULONG;
16 #else
17 typedef long BLASLONG;
18 typedef unsigned long BLASULONG;
19 #endif
20
21 #ifdef LAPACK_ILP64
22 typedef BLASLONG blasint;
23 #if defined(_WIN64)
24 #define blasabs(x) llabs(x)
25 #else
26 #define blasabs(x) labs(x)
27 #endif
28 #else
29 typedef int blasint;
30 #define blasabs(x) abs(x)
31 #endif
32
33 typedef blasint integer;
34
35 typedef unsigned int uinteger;
36 typedef char *address;
37 typedef short int shortint;
38 typedef float real;
39 typedef double doublereal;
40 typedef struct { real r, i; } complex;
41 typedef struct { doublereal r, i; } doublecomplex;
42 #ifdef _MSC_VER
43 static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
44 static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
45 static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
46 static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
47 #else
48 static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
49 static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
50 static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
51 static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
52 #endif
53 #define pCf(z) (*_pCf(z))
54 #define pCd(z) (*_pCd(z))
55 typedef int logical;
56 typedef short int shortlogical;
57 typedef char logical1;
58 typedef char integer1;
59
60 #define TRUE_ (1)
61 #define FALSE_ (0)
62
63 /* Extern is for use with -E */
64 #ifndef Extern
65 #define Extern extern
66 #endif
67
68 /* I/O stuff */
69
70 typedef int flag;
71 typedef int ftnlen;
72 typedef int ftnint;
73
74 /*external read, write*/
75 typedef struct
76 {       flag cierr;
77         ftnint ciunit;
78         flag ciend;
79         char *cifmt;
80         ftnint cirec;
81 } cilist;
82
83 /*internal read, write*/
84 typedef struct
85 {       flag icierr;
86         char *iciunit;
87         flag iciend;
88         char *icifmt;
89         ftnint icirlen;
90         ftnint icirnum;
91 } icilist;
92
93 /*open*/
94 typedef struct
95 {       flag oerr;
96         ftnint ounit;
97         char *ofnm;
98         ftnlen ofnmlen;
99         char *osta;
100         char *oacc;
101         char *ofm;
102         ftnint orl;
103         char *oblnk;
104 } olist;
105
106 /*close*/
107 typedef struct
108 {       flag cerr;
109         ftnint cunit;
110         char *csta;
111 } cllist;
112
113 /*rewind, backspace, endfile*/
114 typedef struct
115 {       flag aerr;
116         ftnint aunit;
117 } alist;
118
119 /* inquire */
120 typedef struct
121 {       flag inerr;
122         ftnint inunit;
123         char *infile;
124         ftnlen infilen;
125         ftnint  *inex;  /*parameters in standard's order*/
126         ftnint  *inopen;
127         ftnint  *innum;
128         ftnint  *innamed;
129         char    *inname;
130         ftnlen  innamlen;
131         char    *inacc;
132         ftnlen  inacclen;
133         char    *inseq;
134         ftnlen  inseqlen;
135         char    *indir;
136         ftnlen  indirlen;
137         char    *infmt;
138         ftnlen  infmtlen;
139         char    *inform;
140         ftnint  informlen;
141         char    *inunf;
142         ftnlen  inunflen;
143         ftnint  *inrecl;
144         ftnint  *innrec;
145         char    *inblank;
146         ftnlen  inblanklen;
147 } inlist;
148
149 #define VOID void
150
151 union Multitype {       /* for multiple entry points */
152         integer1 g;
153         shortint h;
154         integer i;
155         /* longint j; */
156         real r;
157         doublereal d;
158         complex c;
159         doublecomplex z;
160         };
161
162 typedef union Multitype Multitype;
163
164 struct Vardesc {        /* for Namelist */
165         char *name;
166         char *addr;
167         ftnlen *dims;
168         int  type;
169         };
170 typedef struct Vardesc Vardesc;
171
172 struct Namelist {
173         char *name;
174         Vardesc **vars;
175         int nvars;
176         };
177 typedef struct Namelist Namelist;
178
179 #define abs(x) ((x) >= 0 ? (x) : -(x))
180 #define dabs(x) (fabs(x))
181 #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
182 #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
183 #define dmin(a,b) (f2cmin(a,b))
184 #define dmax(a,b) (f2cmax(a,b))
185 #define bit_test(a,b)   ((a) >> (b) & 1)
186 #define bit_clear(a,b)  ((a) & ~((uinteger)1 << (b)))
187 #define bit_set(a,b)    ((a) |  ((uinteger)1 << (b)))
188
189 #define abort_() { sig_die("Fortran abort routine called", 1); }
190 #define c_abs(z) (cabsf(Cf(z)))
191 #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
192 #ifdef _MSC_VER
193 #define c_div(c, a, b) {Cf(c)._Val[0] = (Cf(a)._Val[0]/Cf(b)._Val[0]); Cf(c)._Val[1]=(Cf(a)._Val[1]/Cf(b)._Val[1]);}
194 #define z_div(c, a, b) {Cd(c)._Val[0] = (Cd(a)._Val[0]/Cd(b)._Val[0]); Cd(c)._Val[1]=(Cd(a)._Val[1]/df(b)._Val[1]);}
195 #else
196 #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
197 #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
198 #endif
199 #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
200 #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
201 #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
202 //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
203 #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
204 #define d_abs(x) (fabs(*(x)))
205 #define d_acos(x) (acos(*(x)))
206 #define d_asin(x) (asin(*(x)))
207 #define d_atan(x) (atan(*(x)))
208 #define d_atn2(x, y) (atan2(*(x),*(y)))
209 #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
210 #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
211 #define d_cos(x) (cos(*(x)))
212 #define d_cosh(x) (cosh(*(x)))
213 #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
214 #define d_exp(x) (exp(*(x)))
215 #define d_imag(z) (cimag(Cd(z)))
216 #define r_imag(z) (cimagf(Cf(z)))
217 #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
218 #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
219 #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
220 #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
221 #define d_log(x) (log(*(x)))
222 #define d_mod(x, y) (fmod(*(x), *(y)))
223 #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
224 #define d_nint(x) u_nint(*(x))
225 #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
226 #define d_sign(a,b) u_sign(*(a),*(b))
227 #define r_sign(a,b) u_sign(*(a),*(b))
228 #define d_sin(x) (sin(*(x)))
229 #define d_sinh(x) (sinh(*(x)))
230 #define d_sqrt(x) (sqrt(*(x)))
231 #define d_tan(x) (tan(*(x)))
232 #define d_tanh(x) (tanh(*(x)))
233 #define i_abs(x) abs(*(x))
234 #define i_dnnt(x) ((integer)u_nint(*(x)))
235 #define i_len(s, n) (n)
236 #define i_nint(x) ((integer)u_nint(*(x)))
237 #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
238 #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
239 #define pow_si(B,E) spow_ui(*(B),*(E))
240 #define pow_ri(B,E) spow_ui(*(B),*(E))
241 #define pow_di(B,E) dpow_ui(*(B),*(E))
242 #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
243 #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
244 #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
245 #define s_cat(lpp, rpp, rnp, np, llp) {         ftnlen i, nc, ll; char *f__rp, *lp;     ll = (llp); lp = (lpp);         for(i=0; i < (int)*(np); ++i) {                 nc = ll;                if((rnp)[i] < nc) nc = (rnp)[i];                ll -= nc;               f__rp = (rpp)[i];               while(--nc >= 0) *lp++ = *(f__rp)++;         }  while(--ll >= 0) *lp++ = ' '; }
246 #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
247 #define s_copy(A,B,C,D) { int __i,__m; for (__i=0, __m=f2cmin((C),(D)); __i<__m && (B)[__i] != 0; ++__i) (A)[__i] = (B)[__i]; }
248 #define sig_die(s, kill) { exit(1); }
249 #define s_stop(s, n) {exit(0);}
250 static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
251 #define z_abs(z) (cabs(Cd(z)))
252 #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
253 #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
254 #define myexit_() break;
255 #define mycycle() continue;
256 #define myceiling(w) {ceil(w)}
257 #define myhuge(w) {HUGE_VAL}
258 //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
259 #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
260
261 /* procedure parameter types for -A and -C++ */
262
263 #define F2C_proc_par_types 1
264 #ifdef __cplusplus
265 typedef logical (*L_fp)(...);
266 #else
267 typedef logical (*L_fp)();
268 #endif
269
270 static float spow_ui(float x, integer n) {
271         float pow=1.0; unsigned long int u;
272         if(n != 0) {
273                 if(n < 0) n = -n, x = 1/x;
274                 for(u = n; ; ) {
275                         if(u & 01) pow *= x;
276                         if(u >>= 1) x *= x;
277                         else break;
278                 }
279         }
280         return pow;
281 }
282 static double dpow_ui(double x, integer n) {
283         double pow=1.0; unsigned long int u;
284         if(n != 0) {
285                 if(n < 0) n = -n, x = 1/x;
286                 for(u = n; ; ) {
287                         if(u & 01) pow *= x;
288                         if(u >>= 1) x *= x;
289                         else break;
290                 }
291         }
292         return pow;
293 }
294 #ifdef _MSC_VER
295 static _Fcomplex cpow_ui(complex x, integer n) {
296         complex pow={1.0,0.0}; unsigned long int u;
297                 if(n != 0) {
298                 if(n < 0) n = -n, x.r = 1/x.r, x.i=1/x.i;
299                 for(u = n; ; ) {
300                         if(u & 01) pow.r *= x.r, pow.i *= x.i;
301                         if(u >>= 1) x.r *= x.r, x.i *= x.i;
302                         else break;
303                 }
304         }
305         _Fcomplex p={pow.r, pow.i};
306         return p;
307 }
308 #else
309 static _Complex float cpow_ui(_Complex float x, integer n) {
310         _Complex float pow=1.0; unsigned long int u;
311         if(n != 0) {
312                 if(n < 0) n = -n, x = 1/x;
313                 for(u = n; ; ) {
314                         if(u & 01) pow *= x;
315                         if(u >>= 1) x *= x;
316                         else break;
317                 }
318         }
319         return pow;
320 }
321 #endif
322 #ifdef _MSC_VER
323 static _Dcomplex zpow_ui(_Dcomplex x, integer n) {
324         _Dcomplex pow={1.0,0.0}; unsigned long int u;
325         if(n != 0) {
326                 if(n < 0) n = -n, x._Val[0] = 1/x._Val[0], x._Val[1] =1/x._Val[1];
327                 for(u = n; ; ) {
328                         if(u & 01) pow._Val[0] *= x._Val[0], pow._Val[1] *= x._Val[1];
329                         if(u >>= 1) x._Val[0] *= x._Val[0], x._Val[1] *= x._Val[1];
330                         else break;
331                 }
332         }
333         _Dcomplex p = {pow._Val[0], pow._Val[1]};
334         return p;
335 }
336 #else
337 static _Complex double zpow_ui(_Complex double x, integer n) {
338         _Complex double pow=1.0; unsigned long int u;
339         if(n != 0) {
340                 if(n < 0) n = -n, x = 1/x;
341                 for(u = n; ; ) {
342                         if(u & 01) pow *= x;
343                         if(u >>= 1) x *= x;
344                         else break;
345                 }
346         }
347         return pow;
348 }
349 #endif
350 static integer pow_ii(integer x, integer n) {
351         integer pow; unsigned long int u;
352         if (n <= 0) {
353                 if (n == 0 || x == 1) pow = 1;
354                 else if (x != -1) pow = x == 0 ? 1/x : 0;
355                 else n = -n;
356         }
357         if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
358                 u = n;
359                 for(pow = 1; ; ) {
360                         if(u & 01) pow *= x;
361                         if(u >>= 1) x *= x;
362                         else break;
363                 }
364         }
365         return pow;
366 }
367 static integer dmaxloc_(double *w, integer s, integer e, integer *n)
368 {
369         double m; integer i, mi;
370         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
371                 if (w[i-1]>m) mi=i ,m=w[i-1];
372         return mi-s+1;
373 }
374 static integer smaxloc_(float *w, integer s, integer e, integer *n)
375 {
376         float m; integer i, mi;
377         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
378                 if (w[i-1]>m) mi=i ,m=w[i-1];
379         return mi-s+1;
380 }
381 static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
382         integer n = *n_, incx = *incx_, incy = *incy_, i;
383 #ifdef _MSC_VER
384         _Fcomplex zdotc = {0.0, 0.0};
385         if (incx == 1 && incy == 1) {
386                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
387                         zdotc._Val[0] += conjf(Cf(&x[i]))._Val[0] * Cf(&y[i])._Val[0];
388                         zdotc._Val[1] += conjf(Cf(&x[i]))._Val[1] * Cf(&y[i])._Val[1];
389                 }
390         } else {
391                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
392                         zdotc._Val[0] += conjf(Cf(&x[i*incx]))._Val[0] * Cf(&y[i*incy])._Val[0];
393                         zdotc._Val[1] += conjf(Cf(&x[i*incx]))._Val[1] * Cf(&y[i*incy])._Val[1];
394                 }
395         }
396         pCf(z) = zdotc;
397 }
398 #else
399         _Complex float zdotc = 0.0;
400         if (incx == 1 && incy == 1) {
401                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
402                         zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
403                 }
404         } else {
405                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
406                         zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
407                 }
408         }
409         pCf(z) = zdotc;
410 }
411 #endif
412 static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
413         integer n = *n_, incx = *incx_, incy = *incy_, i;
414 #ifdef _MSC_VER
415         _Dcomplex zdotc = {0.0, 0.0};
416         if (incx == 1 && incy == 1) {
417                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
418                         zdotc._Val[0] += conj(Cd(&x[i]))._Val[0] * Cd(&y[i])._Val[0];
419                         zdotc._Val[1] += conj(Cd(&x[i]))._Val[1] * Cd(&y[i])._Val[1];
420                 }
421         } else {
422                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
423                         zdotc._Val[0] += conj(Cd(&x[i*incx]))._Val[0] * Cd(&y[i*incy])._Val[0];
424                         zdotc._Val[1] += conj(Cd(&x[i*incx]))._Val[1] * Cd(&y[i*incy])._Val[1];
425                 }
426         }
427         pCd(z) = zdotc;
428 }
429 #else
430         _Complex double zdotc = 0.0;
431         if (incx == 1 && incy == 1) {
432                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
433                         zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
434                 }
435         } else {
436                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
437                         zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
438                 }
439         }
440         pCd(z) = zdotc;
441 }
442 #endif  
443 static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
444         integer n = *n_, incx = *incx_, incy = *incy_, i;
445 #ifdef _MSC_VER
446         _Fcomplex zdotc = {0.0, 0.0};
447         if (incx == 1 && incy == 1) {
448                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
449                         zdotc._Val[0] += Cf(&x[i])._Val[0] * Cf(&y[i])._Val[0];
450                         zdotc._Val[1] += Cf(&x[i])._Val[1] * Cf(&y[i])._Val[1];
451                 }
452         } else {
453                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
454                         zdotc._Val[0] += Cf(&x[i*incx])._Val[0] * Cf(&y[i*incy])._Val[0];
455                         zdotc._Val[1] += Cf(&x[i*incx])._Val[1] * Cf(&y[i*incy])._Val[1];
456                 }
457         }
458         pCf(z) = zdotc;
459 }
460 #else
461         _Complex float zdotc = 0.0;
462         if (incx == 1 && incy == 1) {
463                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
464                         zdotc += Cf(&x[i]) * Cf(&y[i]);
465                 }
466         } else {
467                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
468                         zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
469                 }
470         }
471         pCf(z) = zdotc;
472 }
473 #endif
474 static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
475         integer n = *n_, incx = *incx_, incy = *incy_, i;
476 #ifdef _MSC_VER
477         _Dcomplex zdotc = {0.0, 0.0};
478         if (incx == 1 && incy == 1) {
479                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
480                         zdotc._Val[0] += Cd(&x[i])._Val[0] * Cd(&y[i])._Val[0];
481                         zdotc._Val[1] += Cd(&x[i])._Val[1] * Cd(&y[i])._Val[1];
482                 }
483         } else {
484                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
485                         zdotc._Val[0] += Cd(&x[i*incx])._Val[0] * Cd(&y[i*incy])._Val[0];
486                         zdotc._Val[1] += Cd(&x[i*incx])._Val[1] * Cd(&y[i*incy])._Val[1];
487                 }
488         }
489         pCd(z) = zdotc;
490 }
491 #else
492         _Complex double zdotc = 0.0;
493         if (incx == 1 && incy == 1) {
494                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
495                         zdotc += Cd(&x[i]) * Cd(&y[i]);
496                 }
497         } else {
498                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
499                         zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
500                 }
501         }
502         pCd(z) = zdotc;
503 }
504 #endif
505 /*  -- translated by f2c (version 20000121).
506    You must link the resulting object file with the libraries:
507         -lf2c -lm   (in that order)
508 */
509
510
511
512
513 /* Table of constant values */
514
515 static integer c__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> CHEEVR computes the eigenvalues and, optionally, the left and/or right eigenvectors for HE mat
523 rices</b> */
524
525 /*  =========== DOCUMENTATION =========== */
526
527 /* Online html documentation available at */
528 /*            http://www.netlib.org/lapack/explore-html/ */
529
530 /* > \htmlonly */
531 /* > Download CHEEVR + dependencies */
532 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cheevr.
533 f"> */
534 /* > [TGZ]</a> */
535 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cheevr.
536 f"> */
537 /* > [ZIP]</a> */
538 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cheevr.
539 f"> */
540 /* > [TXT]</a> */
541 /* > \endhtmlonly */
542
543 /*  Definition: */
544 /*  =========== */
545
546 /*       SUBROUTINE CHEEVR( JOBZ, RANGE, UPLO, N, A, LDA, VL, VU, IL, IU, */
547 /*                          ABSTOL, M, W, Z, LDZ, ISUPPZ, WORK, LWORK, */
548 /*                          RWORK, LRWORK, IWORK, LIWORK, INFO ) */
549
550 /*       CHARACTER          JOBZ, RANGE, UPLO */
551 /*       INTEGER            IL, INFO, IU, LDA, LDZ, LIWORK, LRWORK, LWORK, */
552 /*      $                   M, N */
553 /*       REAL               ABSTOL, VL, VU */
554 /*       INTEGER            ISUPPZ( * ), IWORK( * ) */
555 /*       REAL               RWORK( * ), W( * ) */
556 /*       COMPLEX            A( LDA, * ), WORK( * ), Z( LDZ, * ) */
557
558
559 /* > \par Purpose: */
560 /*  ============= */
561 /* > */
562 /* > \verbatim */
563 /* > */
564 /* > CHEEVR computes selected eigenvalues and, optionally, eigenvectors */
565 /* > of a complex Hermitian matrix A.  Eigenvalues and eigenvectors can */
566 /* > be selected by specifying either a range of values or a range of */
567 /* > indices for the desired eigenvalues. */
568 /* > */
569 /* > CHEEVR first reduces the matrix A to tridiagonal form T with a call */
570 /* > to CHETRD.  Then, whenever possible, CHEEVR calls CSTEMR to compute */
571 /* > the eigenspectrum using Relatively Robust Representations.  CSTEMR */
572 /* > computes eigenvalues by the dqds algorithm, while orthogonal */
573 /* > eigenvectors are computed from various "good" L D L^T representations */
574 /* > (also known as Relatively Robust Representations). Gram-Schmidt */
575 /* > orthogonalization is avoided as far as possible. More specifically, */
576 /* > the various steps of the algorithm are as follows. */
577 /* > */
578 /* > For each unreduced block (submatrix) of T, */
579 /* >    (a) Compute T - sigma I  = L D L^T, so that L and D */
580 /* >        define all the wanted eigenvalues to high relative accuracy. */
581 /* >        This means that small relative changes in the entries of D and L */
582 /* >        cause only small relative changes in the eigenvalues and */
583 /* >        eigenvectors. The standard (unfactored) representation of the */
584 /* >        tridiagonal matrix T does not have this property in general. */
585 /* >    (b) Compute the eigenvalues to suitable accuracy. */
586 /* >        If the eigenvectors are desired, the algorithm attains full */
587 /* >        accuracy of the computed eigenvalues only right before */
588 /* >        the corresponding vectors have to be computed, see steps c) and d). */
589 /* >    (c) For each cluster of close eigenvalues, select a new */
590 /* >        shift close to the cluster, find a new factorization, and refine */
591 /* >        the shifted eigenvalues to suitable accuracy. */
592 /* >    (d) For each eigenvalue with a large enough relative separation compute */
593 /* >        the corresponding eigenvector by forming a rank revealing twisted */
594 /* >        factorization. Go back to (c) for any clusters that remain. */
595 /* > */
596 /* > The desired accuracy of the output can be specified by the input */
597 /* > parameter ABSTOL. */
598 /* > */
599 /* > For more details, see DSTEMR's documentation and: */
600 /* > - Inderjit S. Dhillon and Beresford N. Parlett: "Multiple representations */
601 /* >   to compute orthogonal eigenvectors of symmetric tridiagonal matrices," */
602 /* >   Linear Algebra and its Applications, 387(1), pp. 1-28, August 2004. */
603 /* > - Inderjit Dhillon and Beresford Parlett: "Orthogonal Eigenvectors and */
604 /* >   Relative Gaps," SIAM Journal on Matrix Analysis and Applications, Vol. 25, */
605 /* >   2004.  Also LAPACK Working Note 154. */
606 /* > - Inderjit Dhillon: "A new O(n^2) algorithm for the symmetric */
607 /* >   tridiagonal eigenvalue/eigenvector problem", */
608 /* >   Computer Science Division Technical Report No. UCB/CSD-97-971, */
609 /* >   UC Berkeley, May 1997. */
610 /* > */
611 /* > */
612 /* > Note 1 : CHEEVR calls CSTEMR when the full spectrum is requested */
613 /* > on machines which conform to the ieee-754 floating point standard. */
614 /* > CHEEVR calls SSTEBZ and CSTEIN on non-ieee machines and */
615 /* > when partial spectrum requests are made. */
616 /* > */
617 /* > Normal execution of CSTEMR may create NaNs and infinities and */
618 /* > hence may abort due to a floating point exception in environments */
619 /* > which do not handle NaNs and infinities in the ieee standard default */
620 /* > manner. */
621 /* > \endverbatim */
622
623 /*  Arguments: */
624 /*  ========== */
625
626 /* > \param[in] JOBZ */
627 /* > \verbatim */
628 /* >          JOBZ is CHARACTER*1 */
629 /* >          = 'N':  Compute eigenvalues only; */
630 /* >          = 'V':  Compute eigenvalues and eigenvectors. */
631 /* > \endverbatim */
632 /* > */
633 /* > \param[in] RANGE */
634 /* > \verbatim */
635 /* >          RANGE is CHARACTER*1 */
636 /* >          = 'A': all eigenvalues will be found. */
637 /* >          = 'V': all eigenvalues in the half-open interval (VL,VU] */
638 /* >                 will be found. */
639 /* >          = 'I': the IL-th through IU-th eigenvalues will be found. */
640 /* >          For RANGE = 'V' or 'I' and IU - IL < N - 1, SSTEBZ and */
641 /* >          CSTEIN are called */
642 /* > \endverbatim */
643 /* > */
644 /* > \param[in] UPLO */
645 /* > \verbatim */
646 /* >          UPLO is CHARACTER*1 */
647 /* >          = 'U':  Upper triangle of A is stored; */
648 /* >          = 'L':  Lower triangle of A is stored. */
649 /* > \endverbatim */
650 /* > */
651 /* > \param[in] N */
652 /* > \verbatim */
653 /* >          N is INTEGER */
654 /* >          The order of the matrix A.  N >= 0. */
655 /* > \endverbatim */
656 /* > */
657 /* > \param[in,out] A */
658 /* > \verbatim */
659 /* >          A is COMPLEX array, dimension (LDA, N) */
660 /* >          On entry, the Hermitian matrix A.  If UPLO = 'U', the */
661 /* >          leading N-by-N upper triangular part of A contains the */
662 /* >          upper triangular part of the matrix A.  If UPLO = 'L', */
663 /* >          the leading N-by-N lower triangular part of A contains */
664 /* >          the lower triangular part of the matrix A. */
665 /* >          On exit, the lower triangle (if UPLO='L') or the upper */
666 /* >          triangle (if UPLO='U') of A, including the diagonal, is */
667 /* >          destroyed. */
668 /* > \endverbatim */
669 /* > */
670 /* > \param[in] LDA */
671 /* > \verbatim */
672 /* >          LDA is INTEGER */
673 /* >          The leading dimension of the array A.  LDA >= f2cmax(1,N). */
674 /* > \endverbatim */
675 /* > */
676 /* > \param[in] VL */
677 /* > \verbatim */
678 /* >          VL is REAL */
679 /* >          If RANGE='V', the lower bound of the interval to */
680 /* >          be searched for eigenvalues. VL < VU. */
681 /* >          Not referenced if RANGE = 'A' or 'I'. */
682 /* > \endverbatim */
683 /* > */
684 /* > \param[in] VU */
685 /* > \verbatim */
686 /* >          VU is REAL */
687 /* >          If RANGE='V', the upper bound of the interval to */
688 /* >          be searched for eigenvalues. VL < VU. */
689 /* >          Not referenced if RANGE = 'A' or 'I'. */
690 /* > \endverbatim */
691 /* > */
692 /* > \param[in] IL */
693 /* > \verbatim */
694 /* >          IL is INTEGER */
695 /* >          If RANGE='I', the index of the */
696 /* >          smallest eigenvalue to be returned. */
697 /* >          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. */
698 /* >          Not referenced if RANGE = 'A' or 'V'. */
699 /* > \endverbatim */
700 /* > */
701 /* > \param[in] IU */
702 /* > \verbatim */
703 /* >          IU is INTEGER */
704 /* >          If RANGE='I', the index of the */
705 /* >          largest eigenvalue to be returned. */
706 /* >          1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0. */
707 /* >          Not referenced if RANGE = 'A' or 'V'. */
708 /* > \endverbatim */
709 /* > */
710 /* > \param[in] ABSTOL */
711 /* > \verbatim */
712 /* >          ABSTOL is REAL */
713 /* >          The absolute error tolerance for the eigenvalues. */
714 /* >          An approximate eigenvalue is accepted as converged */
715 /* >          when it is determined to lie in an interval [a,b] */
716 /* >          of width less than or equal to */
717 /* > */
718 /* >                  ABSTOL + EPS *   f2cmax( |a|,|b| ) , */
719 /* > */
720 /* >          where EPS is the machine precision.  If ABSTOL is less than */
721 /* >          or equal to zero, then  EPS*|T|  will be used in its place, */
722 /* >          where |T| is the 1-norm of the tridiagonal matrix obtained */
723 /* >          by reducing A to tridiagonal form. */
724 /* > */
725 /* >          See "Computing Small Singular Values of Bidiagonal Matrices */
726 /* >          with Guaranteed High Relative Accuracy," by Demmel and */
727 /* >          Kahan, LAPACK Working Note #3. */
728 /* > */
729 /* >          If high relative accuracy is important, set ABSTOL to */
730 /* >          SLAMCH( 'Safe minimum' ).  Doing so will guarantee that */
731 /* >          eigenvalues are computed to high relative accuracy when */
732 /* >          possible in future releases.  The current code does not */
733 /* >          make any guarantees about high relative accuracy, but */
734 /* >          future releases will. See J. Barlow and J. Demmel, */
735 /* >          "Computing Accurate Eigensystems of Scaled Diagonally */
736 /* >          Dominant Matrices", LAPACK Working Note #7, for a discussion */
737 /* >          of which matrices define their eigenvalues to high relative */
738 /* >          accuracy. */
739 /* > \endverbatim */
740 /* > */
741 /* > \param[out] M */
742 /* > \verbatim */
743 /* >          M is INTEGER */
744 /* >          The total number of eigenvalues found.  0 <= M <= N. */
745 /* >          If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1. */
746 /* > \endverbatim */
747 /* > */
748 /* > \param[out] W */
749 /* > \verbatim */
750 /* >          W is REAL array, dimension (N) */
751 /* >          The first M elements contain the selected eigenvalues in */
752 /* >          ascending order. */
753 /* > \endverbatim */
754 /* > */
755 /* > \param[out] Z */
756 /* > \verbatim */
757 /* >          Z is COMPLEX array, dimension (LDZ, f2cmax(1,M)) */
758 /* >          If JOBZ = 'V', then if INFO = 0, the first M columns of Z */
759 /* >          contain the orthonormal eigenvectors of the matrix A */
760 /* >          corresponding to the selected eigenvalues, with the i-th */
761 /* >          column of Z holding the eigenvector associated with W(i). */
762 /* >          If JOBZ = 'N', then Z is not referenced. */
763 /* >          Note: the user must ensure that at least f2cmax(1,M) columns are */
764 /* >          supplied in the array Z; if RANGE = 'V', the exact value of M */
765 /* >          is not known in advance and an upper bound must be used. */
766 /* > \endverbatim */
767 /* > */
768 /* > \param[in] LDZ */
769 /* > \verbatim */
770 /* >          LDZ is INTEGER */
771 /* >          The leading dimension of the array Z.  LDZ >= 1, and if */
772 /* >          JOBZ = 'V', LDZ >= f2cmax(1,N). */
773 /* > \endverbatim */
774 /* > */
775 /* > \param[out] ISUPPZ */
776 /* > \verbatim */
777 /* >          ISUPPZ is INTEGER array, dimension ( 2*f2cmax(1,M) ) */
778 /* >          The support of the eigenvectors in Z, i.e., the indices */
779 /* >          indicating the nonzero elements in Z. The i-th eigenvector */
780 /* >          is nonzero only in elements ISUPPZ( 2*i-1 ) through */
781 /* >          ISUPPZ( 2*i ). This is an output of CSTEMR (tridiagonal */
782 /* >          matrix). The support of the eigenvectors of A is typically */
783 /* >          1:N because of the unitary transformations applied by CUNMTR. */
784 /* >          Implemented only for RANGE = 'A' or 'I' and IU - IL = N - 1 */
785 /* > \endverbatim */
786 /* > */
787 /* > \param[out] WORK */
788 /* > \verbatim */
789 /* >          WORK is COMPLEX array, dimension (MAX(1,LWORK)) */
790 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
791 /* > \endverbatim */
792 /* > */
793 /* > \param[in] LWORK */
794 /* > \verbatim */
795 /* >          LWORK is INTEGER */
796 /* >          The length of the array WORK.  LWORK >= f2cmax(1,2*N). */
797 /* >          For optimal efficiency, LWORK >= (NB+1)*N, */
798 /* >          where NB is the f2cmax of the blocksize for CHETRD and for */
799 /* >          CUNMTR as returned by ILAENV. */
800 /* > */
801 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
802 /* >          only calculates the optimal sizes of the WORK, RWORK and */
803 /* >          IWORK arrays, returns these values as the first entries of */
804 /* >          the WORK, RWORK and IWORK arrays, and no error message */
805 /* >          related to LWORK or LRWORK or LIWORK is issued by XERBLA. */
806 /* > \endverbatim */
807 /* > */
808 /* > \param[out] RWORK */
809 /* > \verbatim */
810 /* >          RWORK is REAL array, dimension (MAX(1,LRWORK)) */
811 /* >          On exit, if INFO = 0, RWORK(1) returns the optimal */
812 /* >          (and minimal) LRWORK. */
813 /* > \endverbatim */
814 /* > */
815 /* > \param[in] LRWORK */
816 /* > \verbatim */
817 /* >          LRWORK is INTEGER */
818 /* >          The length of the array RWORK.  LRWORK >= f2cmax(1,24*N). */
819 /* > */
820 /* >          If LRWORK = -1, then a workspace query is assumed; the */
821 /* >          routine only calculates the optimal sizes of the WORK, RWORK */
822 /* >          and IWORK arrays, returns these values as the first entries */
823 /* >          of the WORK, RWORK and IWORK arrays, and no error message */
824 /* >          related to LWORK or LRWORK or LIWORK is issued by XERBLA. */
825 /* > \endverbatim */
826 /* > */
827 /* > \param[out] IWORK */
828 /* > \verbatim */
829 /* >          IWORK is INTEGER array, dimension (MAX(1,LIWORK)) */
830 /* >          On exit, if INFO = 0, IWORK(1) returns the optimal */
831 /* >          (and minimal) LIWORK. */
832 /* > \endverbatim */
833 /* > */
834 /* > \param[in] LIWORK */
835 /* > \verbatim */
836 /* >          LIWORK is INTEGER */
837 /* >          The dimension of the array IWORK.  LIWORK >= f2cmax(1,10*N). */
838 /* > */
839 /* >          If LIWORK = -1, then a workspace query is assumed; the */
840 /* >          routine only calculates the optimal sizes of the WORK, RWORK */
841 /* >          and IWORK arrays, returns these values as the first entries */
842 /* >          of the WORK, RWORK and IWORK arrays, and no error message */
843 /* >          related to LWORK or LRWORK or LIWORK is issued by XERBLA. */
844 /* > \endverbatim */
845 /* > */
846 /* > \param[out] INFO */
847 /* > \verbatim */
848 /* >          INFO is INTEGER */
849 /* >          = 0:  successful exit */
850 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value */
851 /* >          > 0:  Internal error */
852 /* > \endverbatim */
853
854 /*  Authors: */
855 /*  ======== */
856
857 /* > \author Univ. of Tennessee */
858 /* > \author Univ. of California Berkeley */
859 /* > \author Univ. of Colorado Denver */
860 /* > \author NAG Ltd. */
861
862 /* > \date June 2016 */
863
864 /* > \ingroup complexHEeigen */
865
866 /* > \par Contributors: */
867 /*  ================== */
868 /* > */
869 /* >     Inderjit Dhillon, IBM Almaden, USA \n */
870 /* >     Osni Marques, LBNL/NERSC, USA \n */
871 /* >     Ken Stanley, Computer Science Division, University of */
872 /* >       California at Berkeley, USA \n */
873 /* >     Jason Riedy, Computer Science Division, University of */
874 /* >       California at Berkeley, USA \n */
875 /* > */
876 /*  ===================================================================== */
877 /* Subroutine */ int cheevr_(char *jobz, char *range, char *uplo, integer *n, 
878         complex *a, integer *lda, real *vl, real *vu, integer *il, integer *
879         iu, real *abstol, integer *m, real *w, complex *z__, integer *ldz, 
880         integer *isuppz, complex *work, integer *lwork, real *rwork, integer *
881         lrwork, integer *iwork, integer *liwork, integer *info)
882 {
883     /* System generated locals */
884     integer a_dim1, a_offset, z_dim1, z_offset, i__1, i__2;
885     real r__1, r__2;
886
887     /* Local variables */
888     real anrm;
889     integer imax;
890     real rmin, rmax;
891     logical test;
892     integer itmp1, i__, j, indrd, indre;
893     real sigma;
894     extern logical lsame_(char *, char *);
895     integer iinfo;
896     extern /* Subroutine */ int sscal_(integer *, real *, real *, integer *);
897     char order[1];
898     integer indwk;
899     extern /* Subroutine */ int cswap_(integer *, complex *, integer *, 
900             complex *, integer *);
901     integer lwmin;
902     logical lower;
903     extern /* Subroutine */ int scopy_(integer *, real *, integer *, real *, 
904             integer *);
905     logical wantz;
906     integer nb, jj;
907     logical alleig, indeig;
908     integer iscale, ieeeok, indibl, indrdd, indifl, indree;
909     logical valeig;
910     extern real slamch_(char *);
911     extern /* Subroutine */ int chetrd_(char *, integer *, complex *, integer 
912             *, real *, real *, complex *, complex *, integer *, integer *), csscal_(integer *, real *, complex *, integer *);
913     real safmin;
914     extern integer ilaenv_(integer *, char *, char *, integer *, integer *, 
915             integer *, integer *, ftnlen, ftnlen);
916     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
917     real abstll, bignum;
918     integer indtau, indisp;
919     extern /* Subroutine */ int cstein_(integer *, real *, real *, integer *, 
920             real *, integer *, integer *, complex *, integer *, real *, 
921             integer *, integer *, integer *);
922     integer indiwo, indwkn;
923     extern real clansy_(char *, char *, integer *, complex *, integer *, real 
924             *);
925     extern /* Subroutine */ int cstemr_(char *, char *, integer *, real *, 
926             real *, real *, real *, integer *, integer *, integer *, real *, 
927             complex *, integer *, integer *, integer *, logical *, real *, 
928             integer *, integer *, integer *, integer *);
929     integer indrwk, liwmin;
930     logical tryrac;
931     extern /* Subroutine */ int ssterf_(integer *, real *, real *, integer *);
932     integer lrwmin, llwrkn, llwork, nsplit;
933     real smlnum;
934     extern /* Subroutine */ int cunmtr_(char *, char *, char *, integer *, 
935             integer *, complex *, integer *, complex *, complex *, integer *, 
936             complex *, integer *, integer *), sstebz_(
937             char *, char *, integer *, real *, real *, integer *, integer *, 
938             real *, real *, real *, integer *, integer *, real *, integer *, 
939             integer *, real *, integer *, integer *);
940     logical lquery;
941     integer lwkopt;
942     real eps, vll, vuu;
943     integer llrwork;
944     real tmp1;
945
946
947 /*  -- LAPACK driver routine (version 3.7.0) -- */
948 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
949 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
950 /*     June 2016 */
951
952
953 /* ===================================================================== */
954
955
956 /*     Test the input parameters. */
957
958     /* Parameter adjustments */
959     a_dim1 = *lda;
960     a_offset = 1 + a_dim1 * 1;
961     a -= a_offset;
962     --w;
963     z_dim1 = *ldz;
964     z_offset = 1 + z_dim1 * 1;
965     z__ -= z_offset;
966     --isuppz;
967     --work;
968     --rwork;
969     --iwork;
970
971     /* Function Body */
972     ieeeok = ilaenv_(&c__10, "CHEEVR", "N", &c__1, &c__2, &c__3, &c__4, (
973             ftnlen)6, (ftnlen)1);
974
975     lower = lsame_(uplo, "L");
976     wantz = lsame_(jobz, "V");
977     alleig = lsame_(range, "A");
978     valeig = lsame_(range, "V");
979     indeig = lsame_(range, "I");
980
981     lquery = *lwork == -1 || *lrwork == -1 || *liwork == -1;
982
983 /* Computing MAX */
984     i__1 = 1, i__2 = *n * 24;
985     lrwmin = f2cmax(i__1,i__2);
986 /* Computing MAX */
987     i__1 = 1, i__2 = *n * 10;
988     liwmin = f2cmax(i__1,i__2);
989 /* Computing MAX */
990     i__1 = 1, i__2 = *n << 1;
991     lwmin = f2cmax(i__1,i__2);
992
993     *info = 0;
994     if (! (wantz || lsame_(jobz, "N"))) {
995         *info = -1;
996     } else if (! (alleig || valeig || indeig)) {
997         *info = -2;
998     } else if (! (lower || lsame_(uplo, "U"))) {
999         *info = -3;
1000     } else if (*n < 0) {
1001         *info = -4;
1002     } else if (*lda < f2cmax(1,*n)) {
1003         *info = -6;
1004     } else {
1005         if (valeig) {
1006             if (*n > 0 && *vu <= *vl) {
1007                 *info = -8;
1008             }
1009         } else if (indeig) {
1010             if (*il < 1 || *il > f2cmax(1,*n)) {
1011                 *info = -9;
1012             } else if (*iu < f2cmin(*n,*il) || *iu > *n) {
1013                 *info = -10;
1014             }
1015         }
1016     }
1017     if (*info == 0) {
1018         if (*ldz < 1 || wantz && *ldz < *n) {
1019             *info = -15;
1020         }
1021     }
1022
1023     if (*info == 0) {
1024         nb = ilaenv_(&c__1, "CHETRD", uplo, n, &c_n1, &c_n1, &c_n1, (ftnlen)6,
1025                  (ftnlen)1);
1026 /* Computing MAX */
1027         i__1 = nb, i__2 = ilaenv_(&c__1, "CUNMTR", uplo, n, &c_n1, &c_n1, &
1028                 c_n1, (ftnlen)6, (ftnlen)1);
1029         nb = f2cmax(i__1,i__2);
1030 /* Computing MAX */
1031         i__1 = (nb + 1) * *n;
1032         lwkopt = f2cmax(i__1,lwmin);
1033         work[1].r = (real) lwkopt, work[1].i = 0.f;
1034         rwork[1] = (real) lrwmin;
1035         iwork[1] = liwmin;
1036
1037         if (*lwork < lwmin && ! lquery) {
1038             *info = -18;
1039         } else if (*lrwork < lrwmin && ! lquery) {
1040             *info = -20;
1041         } else if (*liwork < liwmin && ! lquery) {
1042             *info = -22;
1043         }
1044     }
1045
1046     if (*info != 0) {
1047         i__1 = -(*info);
1048         xerbla_("CHEEVR", &i__1, (ftnlen)6);
1049         return 0;
1050     } else if (lquery) {
1051         return 0;
1052     }
1053
1054 /*     Quick return if possible */
1055
1056     *m = 0;
1057     if (*n == 0) {
1058         work[1].r = 1.f, work[1].i = 0.f;
1059         return 0;
1060     }
1061
1062     if (*n == 1) {
1063         work[1].r = 2.f, work[1].i = 0.f;
1064         if (alleig || indeig) {
1065             *m = 1;
1066             i__1 = a_dim1 + 1;
1067             w[1] = a[i__1].r;
1068         } else {
1069             i__1 = a_dim1 + 1;
1070             i__2 = a_dim1 + 1;
1071             if (*vl < a[i__1].r && *vu >= a[i__2].r) {
1072                 *m = 1;
1073                 i__1 = a_dim1 + 1;
1074                 w[1] = a[i__1].r;
1075             }
1076         }
1077         if (wantz) {
1078             i__1 = z_dim1 + 1;
1079             z__[i__1].r = 1.f, z__[i__1].i = 0.f;
1080             isuppz[1] = 1;
1081             isuppz[2] = 1;
1082         }
1083         return 0;
1084     }
1085
1086 /*     Get machine constants. */
1087
1088     safmin = slamch_("Safe minimum");
1089     eps = slamch_("Precision");
1090     smlnum = safmin / eps;
1091     bignum = 1.f / smlnum;
1092     rmin = sqrt(smlnum);
1093 /* Computing MIN */
1094     r__1 = sqrt(bignum), r__2 = 1.f / sqrt(sqrt(safmin));
1095     rmax = f2cmin(r__1,r__2);
1096
1097 /*     Scale matrix to allowable range, if necessary. */
1098
1099     iscale = 0;
1100     abstll = *abstol;
1101     if (valeig) {
1102         vll = *vl;
1103         vuu = *vu;
1104     }
1105     anrm = clansy_("M", uplo, n, &a[a_offset], lda, &rwork[1]);
1106     if (anrm > 0.f && anrm < rmin) {
1107         iscale = 1;
1108         sigma = rmin / anrm;
1109     } else if (anrm > rmax) {
1110         iscale = 1;
1111         sigma = rmax / anrm;
1112     }
1113     if (iscale == 1) {
1114         if (lower) {
1115             i__1 = *n;
1116             for (j = 1; j <= i__1; ++j) {
1117                 i__2 = *n - j + 1;
1118                 csscal_(&i__2, &sigma, &a[j + j * a_dim1], &c__1);
1119 /* L10: */
1120             }
1121         } else {
1122             i__1 = *n;
1123             for (j = 1; j <= i__1; ++j) {
1124                 csscal_(&j, &sigma, &a[j * a_dim1 + 1], &c__1);
1125 /* L20: */
1126             }
1127         }
1128         if (*abstol > 0.f) {
1129             abstll = *abstol * sigma;
1130         }
1131         if (valeig) {
1132             vll = *vl * sigma;
1133             vuu = *vu * sigma;
1134         }
1135     }
1136 /*     Initialize indices into workspaces.  Note: The IWORK indices are */
1137 /*     used only if SSTERF or CSTEMR fail. */
1138 /*     WORK(INDTAU:INDTAU+N-1) stores the complex scalar factors of the */
1139 /*     elementary reflectors used in CHETRD. */
1140     indtau = 1;
1141 /*     INDWK is the starting offset of the remaining complex workspace, */
1142 /*     and LLWORK is the remaining complex workspace size. */
1143     indwk = indtau + *n;
1144     llwork = *lwork - indwk + 1;
1145 /*     RWORK(INDRD:INDRD+N-1) stores the real tridiagonal's diagonal */
1146 /*     entries. */
1147     indrd = 1;
1148 /*     RWORK(INDRE:INDRE+N-1) stores the off-diagonal entries of the */
1149 /*     tridiagonal matrix from CHETRD. */
1150     indre = indrd + *n;
1151 /*     RWORK(INDRDD:INDRDD+N-1) is a copy of the diagonal entries over */
1152 /*     -written by CSTEMR (the SSTERF path copies the diagonal to W). */
1153     indrdd = indre + *n;
1154 /*     RWORK(INDREE:INDREE+N-1) is a copy of the off-diagonal entries over */
1155 /*     -written while computing the eigenvalues in SSTERF and CSTEMR. */
1156     indree = indrdd + *n;
1157 /*     INDRWK is the starting offset of the left-over real workspace, and */
1158 /*     LLRWORK is the remaining workspace size. */
1159     indrwk = indree + *n;
1160     llrwork = *lrwork - indrwk + 1;
1161 /*     IWORK(INDIBL:INDIBL+M-1) corresponds to IBLOCK in SSTEBZ and */
1162 /*     stores the block indices of each of the M<=N eigenvalues. */
1163     indibl = 1;
1164 /*     IWORK(INDISP:INDISP+NSPLIT-1) corresponds to ISPLIT in SSTEBZ and */
1165 /*     stores the starting and finishing indices of each block. */
1166     indisp = indibl + *n;
1167 /*     IWORK(INDIFL:INDIFL+N-1) stores the indices of eigenvectors */
1168 /*     that corresponding to eigenvectors that fail to converge in */
1169 /*     SSTEIN.  This information is discarded; if any fail, the driver */
1170 /*     returns INFO > 0. */
1171     indifl = indisp + *n;
1172 /*     INDIWO is the offset of the remaining integer workspace. */
1173     indiwo = indifl + *n;
1174
1175 /*     Call CHETRD to reduce Hermitian matrix to tridiagonal form. */
1176
1177     chetrd_(uplo, n, &a[a_offset], lda, &rwork[indrd], &rwork[indre], &work[
1178             indtau], &work[indwk], &llwork, &iinfo);
1179
1180 /*     If all eigenvalues are desired */
1181 /*     then call SSTERF or CSTEMR and CUNMTR. */
1182
1183     test = FALSE_;
1184     if (indeig) {
1185         if (*il == 1 && *iu == *n) {
1186             test = TRUE_;
1187         }
1188     }
1189     if ((alleig || test) && ieeeok == 1) {
1190         if (! wantz) {
1191             scopy_(n, &rwork[indrd], &c__1, &w[1], &c__1);
1192             i__1 = *n - 1;
1193             scopy_(&i__1, &rwork[indre], &c__1, &rwork[indree], &c__1);
1194             ssterf_(n, &w[1], &rwork[indree], info);
1195         } else {
1196             i__1 = *n - 1;
1197             scopy_(&i__1, &rwork[indre], &c__1, &rwork[indree], &c__1);
1198             scopy_(n, &rwork[indrd], &c__1, &rwork[indrdd], &c__1);
1199
1200             if (*abstol <= *n * 2.f * eps) {
1201                 tryrac = TRUE_;
1202             } else {
1203                 tryrac = FALSE_;
1204             }
1205             cstemr_(jobz, "A", n, &rwork[indrdd], &rwork[indree], vl, vu, il, 
1206                     iu, m, &w[1], &z__[z_offset], ldz, n, &isuppz[1], &tryrac,
1207                      &rwork[indrwk], &llrwork, &iwork[1], liwork, info);
1208
1209 /*           Apply unitary matrix used in reduction to tridiagonal */
1210 /*           form to eigenvectors returned by CSTEMR. */
1211
1212             if (wantz && *info == 0) {
1213                 indwkn = indwk;
1214                 llwrkn = *lwork - indwkn + 1;
1215                 cunmtr_("L", uplo, "N", n, m, &a[a_offset], lda, &work[indtau]
1216                         , &z__[z_offset], ldz, &work[indwkn], &llwrkn, &iinfo);
1217             }
1218         }
1219
1220
1221         if (*info == 0) {
1222             *m = *n;
1223             goto L30;
1224         }
1225         *info = 0;
1226     }
1227
1228 /*     Otherwise, call SSTEBZ and, if eigenvectors are desired, CSTEIN. */
1229 /*     Also call SSTEBZ and CSTEIN if CSTEMR fails. */
1230
1231     if (wantz) {
1232         *(unsigned char *)order = 'B';
1233     } else {
1234         *(unsigned char *)order = 'E';
1235     }
1236     sstebz_(range, order, n, &vll, &vuu, il, iu, &abstll, &rwork[indrd], &
1237             rwork[indre], m, &nsplit, &w[1], &iwork[indibl], &iwork[indisp], &
1238             rwork[indrwk], &iwork[indiwo], info);
1239
1240     if (wantz) {
1241         cstein_(n, &rwork[indrd], &rwork[indre], m, &w[1], &iwork[indibl], &
1242                 iwork[indisp], &z__[z_offset], ldz, &rwork[indrwk], &iwork[
1243                 indiwo], &iwork[indifl], info);
1244
1245 /*        Apply unitary matrix used in reduction to tridiagonal */
1246 /*        form to eigenvectors returned by CSTEIN. */
1247
1248         indwkn = indwk;
1249         llwrkn = *lwork - indwkn + 1;
1250         cunmtr_("L", uplo, "N", n, m, &a[a_offset], lda, &work[indtau], &z__[
1251                 z_offset], ldz, &work[indwkn], &llwrkn, &iinfo);
1252     }
1253
1254 /*     If matrix was scaled, then rescale eigenvalues appropriately. */
1255
1256 L30:
1257     if (iscale == 1) {
1258         if (*info == 0) {
1259             imax = *m;
1260         } else {
1261             imax = *info - 1;
1262         }
1263         r__1 = 1.f / sigma;
1264         sscal_(&imax, &r__1, &w[1], &c__1);
1265     }
1266
1267 /*     If eigenvalues are not in order, then sort them, along with */
1268 /*     eigenvectors. */
1269
1270     if (wantz) {
1271         i__1 = *m - 1;
1272         for (j = 1; j <= i__1; ++j) {
1273             i__ = 0;
1274             tmp1 = w[j];
1275             i__2 = *m;
1276             for (jj = j + 1; jj <= i__2; ++jj) {
1277                 if (w[jj] < tmp1) {
1278                     i__ = jj;
1279                     tmp1 = w[jj];
1280                 }
1281 /* L40: */
1282             }
1283
1284             if (i__ != 0) {
1285                 itmp1 = iwork[indibl + i__ - 1];
1286                 w[i__] = w[j];
1287                 iwork[indibl + i__ - 1] = iwork[indibl + j - 1];
1288                 w[j] = tmp1;
1289                 iwork[indibl + j - 1] = itmp1;
1290                 cswap_(n, &z__[i__ * z_dim1 + 1], &c__1, &z__[j * z_dim1 + 1],
1291                          &c__1);
1292             }
1293 /* L50: */
1294         }
1295     }
1296
1297 /*     Set WORK(1) to optimal workspace size. */
1298
1299     work[1].r = (real) lwkopt, work[1].i = 0.f;
1300     rwork[1] = (real) lrwmin;
1301     iwork[1] = liwmin;
1302
1303     return 0;
1304
1305 /*     End of CHEEVR */
1306
1307 } /* cheevr_ */
1308