C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / dlarre.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__1 = 1;
516 static integer c__2 = 2;
517
518 /* > \brief \b DLARRE given the tridiagonal matrix T, sets small off-diagonal elements to zero and for each un
519 reduced block Ti, finds base representations and eigenvalues. */
520
521 /*  =========== DOCUMENTATION =========== */
522
523 /* Online html documentation available at */
524 /*            http://www.netlib.org/lapack/explore-html/ */
525
526 /* > \htmlonly */
527 /* > Download DLARRE + dependencies */
528 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlarre.
529 f"> */
530 /* > [TGZ]</a> */
531 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlarre.
532 f"> */
533 /* > [ZIP]</a> */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlarre.
535 f"> */
536 /* > [TXT]</a> */
537 /* > \endhtmlonly */
538
539 /*  Definition: */
540 /*  =========== */
541
542 /*       SUBROUTINE DLARRE( RANGE, N, VL, VU, IL, IU, D, E, E2, */
543 /*                           RTOL1, RTOL2, SPLTOL, NSPLIT, ISPLIT, M, */
544 /*                           W, WERR, WGAP, IBLOCK, INDEXW, GERS, PIVMIN, */
545 /*                           WORK, IWORK, INFO ) */
546
547 /*       CHARACTER          RANGE */
548 /*       INTEGER            IL, INFO, IU, M, N, NSPLIT */
549 /*       DOUBLE PRECISION  PIVMIN, RTOL1, RTOL2, SPLTOL, VL, VU */
550 /*       INTEGER            IBLOCK( * ), ISPLIT( * ), IWORK( * ), */
551 /*      $                   INDEXW( * ) */
552 /*       DOUBLE PRECISION   D( * ), E( * ), E2( * ), GERS( * ), */
553 /*      $                   W( * ),WERR( * ), WGAP( * ), WORK( * ) */
554
555
556 /* > \par Purpose: */
557 /*  ============= */
558 /* > */
559 /* > \verbatim */
560 /* > */
561 /* > To find the desired eigenvalues of a given real symmetric */
562 /* > tridiagonal matrix T, DLARRE sets any "small" off-diagonal */
563 /* > elements to zero, and for each unreduced block T_i, it finds */
564 /* > (a) a suitable shift at one end of the block's spectrum, */
565 /* > (b) the base representation, T_i - sigma_i I = L_i D_i L_i^T, and */
566 /* > (c) eigenvalues of each L_i D_i L_i^T. */
567 /* > The representations and eigenvalues found are then used by */
568 /* > DSTEMR to compute the eigenvectors of T. */
569 /* > The accuracy varies depending on whether bisection is used to */
570 /* > find a few eigenvalues or the dqds algorithm (subroutine DLASQ2) to */
571 /* > conpute all and then discard any unwanted one. */
572 /* > As an added benefit, DLARRE also outputs the n */
573 /* > Gerschgorin intervals for the matrices L_i D_i L_i^T. */
574 /* > \endverbatim */
575
576 /*  Arguments: */
577 /*  ========== */
578
579 /* > \param[in] RANGE */
580 /* > \verbatim */
581 /* >          RANGE is CHARACTER*1 */
582 /* >          = 'A': ("All")   all eigenvalues will be found. */
583 /* >          = 'V': ("Value") all eigenvalues in the half-open interval */
584 /* >                           (VL, VU] will be found. */
585 /* >          = 'I': ("Index") the IL-th through IU-th eigenvalues (of the */
586 /* >                           entire matrix) will be found. */
587 /* > \endverbatim */
588 /* > */
589 /* > \param[in] N */
590 /* > \verbatim */
591 /* >          N is INTEGER */
592 /* >          The order of the matrix. N > 0. */
593 /* > \endverbatim */
594 /* > */
595 /* > \param[in,out] VL */
596 /* > \verbatim */
597 /* >          VL is DOUBLE PRECISION */
598 /* >          If RANGE='V', the lower bound for the eigenvalues. */
599 /* >          Eigenvalues less than or equal to VL, or greater than VU, */
600 /* >          will not be returned.  VL < VU. */
601 /* >          If RANGE='I' or ='A', DLARRE computes bounds on the desired */
602 /* >          part of the spectrum. */
603 /* > \endverbatim */
604 /* > */
605 /* > \param[in,out] VU */
606 /* > \verbatim */
607 /* >          VU is DOUBLE PRECISION */
608 /* >          If RANGE='V', the upper bound for the eigenvalues. */
609 /* >          Eigenvalues less than or equal to VL, or greater than VU, */
610 /* >          will not be returned.  VL < VU. */
611 /* >          If RANGE='I' or ='A', DLARRE computes bounds on the desired */
612 /* >          part of the spectrum. */
613 /* > \endverbatim */
614 /* > */
615 /* > \param[in] IL */
616 /* > \verbatim */
617 /* >          IL is INTEGER */
618 /* >          If RANGE='I', the index of the */
619 /* >          smallest eigenvalue to be returned. */
620 /* >          1 <= IL <= IU <= N. */
621 /* > \endverbatim */
622 /* > */
623 /* > \param[in] IU */
624 /* > \verbatim */
625 /* >          IU is INTEGER */
626 /* >          If RANGE='I', the index of the */
627 /* >          largest eigenvalue to be returned. */
628 /* >          1 <= IL <= IU <= N. */
629 /* > \endverbatim */
630 /* > */
631 /* > \param[in,out] D */
632 /* > \verbatim */
633 /* >          D is DOUBLE PRECISION array, dimension (N) */
634 /* >          On entry, the N diagonal elements of the tridiagonal */
635 /* >          matrix T. */
636 /* >          On exit, the N diagonal elements of the diagonal */
637 /* >          matrices D_i. */
638 /* > \endverbatim */
639 /* > */
640 /* > \param[in,out] E */
641 /* > \verbatim */
642 /* >          E is DOUBLE PRECISION array, dimension (N) */
643 /* >          On entry, the first (N-1) entries contain the subdiagonal */
644 /* >          elements of the tridiagonal matrix T; E(N) need not be set. */
645 /* >          On exit, E contains the subdiagonal elements of the unit */
646 /* >          bidiagonal matrices L_i. The entries E( ISPLIT( I ) ), */
647 /* >          1 <= I <= NSPLIT, contain the base points sigma_i on output. */
648 /* > \endverbatim */
649 /* > */
650 /* > \param[in,out] E2 */
651 /* > \verbatim */
652 /* >          E2 is DOUBLE PRECISION array, dimension (N) */
653 /* >          On entry, the first (N-1) entries contain the SQUARES of the */
654 /* >          subdiagonal elements of the tridiagonal matrix T; */
655 /* >          E2(N) need not be set. */
656 /* >          On exit, the entries E2( ISPLIT( I ) ), */
657 /* >          1 <= I <= NSPLIT, have been set to zero */
658 /* > \endverbatim */
659 /* > */
660 /* > \param[in] RTOL1 */
661 /* > \verbatim */
662 /* >          RTOL1 is DOUBLE PRECISION */
663 /* > \endverbatim */
664 /* > */
665 /* > \param[in] RTOL2 */
666 /* > \verbatim */
667 /* >          RTOL2 is DOUBLE PRECISION */
668 /* >           Parameters for bisection. */
669 /* >           An interval [LEFT,RIGHT] has converged if */
670 /* >           RIGHT-LEFT < MAX( RTOL1*GAP, RTOL2*MAX(|LEFT|,|RIGHT|) ) */
671 /* > \endverbatim */
672 /* > */
673 /* > \param[in] SPLTOL */
674 /* > \verbatim */
675 /* >          SPLTOL is DOUBLE PRECISION */
676 /* >          The threshold for splitting. */
677 /* > \endverbatim */
678 /* > */
679 /* > \param[out] NSPLIT */
680 /* > \verbatim */
681 /* >          NSPLIT is INTEGER */
682 /* >          The number of blocks T splits into. 1 <= NSPLIT <= N. */
683 /* > \endverbatim */
684 /* > */
685 /* > \param[out] ISPLIT */
686 /* > \verbatim */
687 /* >          ISPLIT is INTEGER array, dimension (N) */
688 /* >          The splitting points, at which T breaks up into blocks. */
689 /* >          The first block consists of rows/columns 1 to ISPLIT(1), */
690 /* >          the second of rows/columns ISPLIT(1)+1 through ISPLIT(2), */
691 /* >          etc., and the NSPLIT-th consists of rows/columns */
692 /* >          ISPLIT(NSPLIT-1)+1 through ISPLIT(NSPLIT)=N. */
693 /* > \endverbatim */
694 /* > */
695 /* > \param[out] M */
696 /* > \verbatim */
697 /* >          M is INTEGER */
698 /* >          The total number of eigenvalues (of all L_i D_i L_i^T) */
699 /* >          found. */
700 /* > \endverbatim */
701 /* > */
702 /* > \param[out] W */
703 /* > \verbatim */
704 /* >          W is DOUBLE PRECISION array, dimension (N) */
705 /* >          The first M elements contain the eigenvalues. The */
706 /* >          eigenvalues of each of the blocks, L_i D_i L_i^T, are */
707 /* >          sorted in ascending order ( DLARRE may use the */
708 /* >          remaining N-M elements as workspace). */
709 /* > \endverbatim */
710 /* > */
711 /* > \param[out] WERR */
712 /* > \verbatim */
713 /* >          WERR is DOUBLE PRECISION array, dimension (N) */
714 /* >          The error bound on the corresponding eigenvalue in W. */
715 /* > \endverbatim */
716 /* > */
717 /* > \param[out] WGAP */
718 /* > \verbatim */
719 /* >          WGAP is DOUBLE PRECISION array, dimension (N) */
720 /* >          The separation from the right neighbor eigenvalue in W. */
721 /* >          The gap is only with respect to the eigenvalues of the same block */
722 /* >          as each block has its own representation tree. */
723 /* >          Exception: at the right end of a block we store the left gap */
724 /* > \endverbatim */
725 /* > */
726 /* > \param[out] IBLOCK */
727 /* > \verbatim */
728 /* >          IBLOCK is INTEGER array, dimension (N) */
729 /* >          The indices of the blocks (submatrices) associated with the */
730 /* >          corresponding eigenvalues in W; IBLOCK(i)=1 if eigenvalue */
731 /* >          W(i) belongs to the first block from the top, =2 if W(i) */
732 /* >          belongs to the second block, etc. */
733 /* > \endverbatim */
734 /* > */
735 /* > \param[out] INDEXW */
736 /* > \verbatim */
737 /* >          INDEXW is INTEGER array, dimension (N) */
738 /* >          The indices of the eigenvalues within each block (submatrix); */
739 /* >          for example, INDEXW(i)= 10 and IBLOCK(i)=2 imply that the */
740 /* >          i-th eigenvalue W(i) is the 10-th eigenvalue in block 2 */
741 /* > \endverbatim */
742 /* > */
743 /* > \param[out] GERS */
744 /* > \verbatim */
745 /* >          GERS is DOUBLE PRECISION array, dimension (2*N) */
746 /* >          The N Gerschgorin intervals (the i-th Gerschgorin interval */
747 /* >          is (GERS(2*i-1), GERS(2*i)). */
748 /* > \endverbatim */
749 /* > */
750 /* > \param[out] PIVMIN */
751 /* > \verbatim */
752 /* >          PIVMIN is DOUBLE PRECISION */
753 /* >          The minimum pivot in the Sturm sequence for T. */
754 /* > \endverbatim */
755 /* > */
756 /* > \param[out] WORK */
757 /* > \verbatim */
758 /* >          WORK is DOUBLE PRECISION array, dimension (6*N) */
759 /* >          Workspace. */
760 /* > \endverbatim */
761 /* > */
762 /* > \param[out] IWORK */
763 /* > \verbatim */
764 /* >          IWORK is INTEGER array, dimension (5*N) */
765 /* >          Workspace. */
766 /* > \endverbatim */
767 /* > */
768 /* > \param[out] INFO */
769 /* > \verbatim */
770 /* >          INFO is INTEGER */
771 /* >          = 0:  successful exit */
772 /* >          > 0:  A problem occurred in DLARRE. */
773 /* >          < 0:  One of the called subroutines signaled an internal problem. */
774 /* >                Needs inspection of the corresponding parameter IINFO */
775 /* >                for further information. */
776 /* > */
777 /* >          =-1:  Problem in DLARRD. */
778 /* >          = 2:  No base representation could be found in MAXTRY iterations. */
779 /* >                Increasing MAXTRY and recompilation might be a remedy. */
780 /* >          =-3:  Problem in DLARRB when computing the refined root */
781 /* >                representation for DLASQ2. */
782 /* >          =-4:  Problem in DLARRB when preforming bisection on the */
783 /* >                desired part of the spectrum. */
784 /* >          =-5:  Problem in DLASQ2. */
785 /* >          =-6:  Problem in DLASQ2. */
786 /* > \endverbatim */
787
788 /*  Authors: */
789 /*  ======== */
790
791 /* > \author Univ. of Tennessee */
792 /* > \author Univ. of California Berkeley */
793 /* > \author Univ. of Colorado Denver */
794 /* > \author NAG Ltd. */
795
796 /* > \date June 2016 */
797
798 /* > \ingroup OTHERauxiliary */
799
800 /* > \par Further Details: */
801 /*  ===================== */
802 /* > */
803 /* > \verbatim */
804 /* > */
805 /* >  The base representations are required to suffer very little */
806 /* >  element growth and consequently define all their eigenvalues to */
807 /* >  high relative accuracy. */
808 /* > \endverbatim */
809
810 /* > \par Contributors: */
811 /*  ================== */
812 /* > */
813 /* >     Beresford Parlett, University of California, Berkeley, USA \n */
814 /* >     Jim Demmel, University of California, Berkeley, USA \n */
815 /* >     Inderjit Dhillon, University of Texas, Austin, USA \n */
816 /* >     Osni Marques, LBNL/NERSC, USA \n */
817 /* >     Christof Voemel, University of California, Berkeley, USA \n */
818 /* > */
819 /*  ===================================================================== */
820 /* Subroutine */ int dlarre_(char *range, integer *n, doublereal *vl, 
821         doublereal *vu, integer *il, integer *iu, doublereal *d__, doublereal 
822         *e, doublereal *e2, doublereal *rtol1, doublereal *rtol2, doublereal *
823         spltol, integer *nsplit, integer *isplit, integer *m, doublereal *w, 
824         doublereal *werr, doublereal *wgap, integer *iblock, integer *indexw, 
825         doublereal *gers, doublereal *pivmin, doublereal *work, integer *
826         iwork, integer *info)
827 {
828     /* System generated locals */
829     integer i__1, i__2;
830     doublereal d__1, d__2, d__3;
831
832     /* Local variables */
833     doublereal eabs;
834     integer iend, jblk;
835     doublereal eold;
836     integer indl;
837     doublereal dmax__, emax;
838     integer wend, idum, indu;
839     doublereal rtol;
840     integer i__, j, iseed[4];
841     doublereal avgap, sigma;
842     extern logical lsame_(char *, char *);
843     integer iinfo;
844     extern /* Subroutine */ int dcopy_(integer *, doublereal *, integer *, 
845             doublereal *, integer *);
846     logical norep;
847     doublereal s1, s2;
848     extern /* Subroutine */ int dlasq2_(integer *, doublereal *, integer *);
849     integer mb;
850     doublereal gl;
851     integer in;
852     extern doublereal dlamch_(char *);
853     integer mm;
854     doublereal gu;
855     integer ibegin;
856     logical forceb;
857     integer irange;
858     doublereal sgndef;
859     extern /* Subroutine */ int dlarra_(integer *, doublereal *, doublereal *,
860              doublereal *, doublereal *, doublereal *, integer *, integer *, 
861             integer *), dlarrb_(integer *, doublereal *, doublereal *, 
862             integer *, integer *, doublereal *, doublereal *, integer *, 
863             doublereal *, doublereal *, doublereal *, doublereal *, integer *,
864              doublereal *, doublereal *, integer *, integer *), dlarrc_(char *
865             , integer *, doublereal *, doublereal *, doublereal *, doublereal 
866             *, doublereal *, integer *, integer *, integer *, integer *);
867     integer wbegin;
868     doublereal safmin, spdiam;
869     extern /* Subroutine */ int dlarrd_(char *, char *, integer *, doublereal 
870             *, doublereal *, integer *, integer *, doublereal *, doublereal *,
871              doublereal *, doublereal *, doublereal *, doublereal *, integer *
872             , integer *, integer *, doublereal *, doublereal *, doublereal *, 
873             doublereal *, integer *, integer *, doublereal *, integer *, 
874             integer *), dlarrk_(integer *, integer *, 
875             doublereal *, doublereal *, doublereal *, doublereal *, 
876             doublereal *, doublereal *, doublereal *, doublereal *, integer *)
877             ;
878     logical usedqd;
879     doublereal clwdth, isleft;
880     extern /* Subroutine */ int dlarnv_(integer *, integer *, integer *, 
881             doublereal *);
882     doublereal isrght, bsrtol, dpivot;
883     integer cnt;
884     doublereal eps, tau, tmp, rtl;
885     integer cnt1, cnt2;
886     doublereal tmp1;
887
888
889 /*  -- LAPACK auxiliary routine (version 3.8.0) -- */
890 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
891 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
892 /*     June 2016 */
893
894
895 /*  ===================================================================== */
896
897
898     /* Parameter adjustments */
899     --iwork;
900     --work;
901     --gers;
902     --indexw;
903     --iblock;
904     --wgap;
905     --werr;
906     --w;
907     --isplit;
908     --e2;
909     --e;
910     --d__;
911
912     /* Function Body */
913     *info = 0;
914
915 /*     Quick return if possible */
916
917     if (*n <= 0) {
918         return 0;
919     }
920
921 /*     Decode RANGE */
922
923     if (lsame_(range, "A")) {
924         irange = 1;
925     } else if (lsame_(range, "V")) {
926         irange = 3;
927     } else if (lsame_(range, "I")) {
928         irange = 2;
929     }
930     *m = 0;
931 /*     Get machine constants */
932     safmin = dlamch_("S");
933     eps = dlamch_("P");
934 /*     Set parameters */
935     rtl = sqrt(eps);
936     bsrtol = sqrt(eps);
937 /*     Treat case of 1x1 matrix for quick return */
938     if (*n == 1) {
939         if (irange == 1 || irange == 3 && d__[1] > *vl && d__[1] <= *vu || 
940                 irange == 2 && *il == 1 && *iu == 1) {
941             *m = 1;
942             w[1] = d__[1];
943 /*           The computation error of the eigenvalue is zero */
944             werr[1] = 0.;
945             wgap[1] = 0.;
946             iblock[1] = 1;
947             indexw[1] = 1;
948             gers[1] = d__[1];
949             gers[2] = d__[1];
950         }
951 /*        store the shift for the initial RRR, which is zero in this case */
952         e[1] = 0.;
953         return 0;
954     }
955 /*     General case: tridiagonal matrix of order > 1 */
956
957 /*     Init WERR, WGAP. Compute Gerschgorin intervals and spectral diameter. */
958 /*     Compute maximum off-diagonal entry and pivmin. */
959     gl = d__[1];
960     gu = d__[1];
961     eold = 0.;
962     emax = 0.;
963     e[*n] = 0.;
964     i__1 = *n;
965     for (i__ = 1; i__ <= i__1; ++i__) {
966         werr[i__] = 0.;
967         wgap[i__] = 0.;
968         eabs = (d__1 = e[i__], abs(d__1));
969         if (eabs >= emax) {
970             emax = eabs;
971         }
972         tmp1 = eabs + eold;
973         gers[(i__ << 1) - 1] = d__[i__] - tmp1;
974 /* Computing MIN */
975         d__1 = gl, d__2 = gers[(i__ << 1) - 1];
976         gl = f2cmin(d__1,d__2);
977         gers[i__ * 2] = d__[i__] + tmp1;
978 /* Computing MAX */
979         d__1 = gu, d__2 = gers[i__ * 2];
980         gu = f2cmax(d__1,d__2);
981         eold = eabs;
982 /* L5: */
983     }
984 /*     The minimum pivot allowed in the Sturm sequence for T */
985 /* Computing MAX */
986 /* Computing 2nd power */
987     d__3 = emax;
988     d__1 = 1., d__2 = d__3 * d__3;
989     *pivmin = safmin * f2cmax(d__1,d__2);
990 /*     Compute spectral diameter. The Gerschgorin bounds give an */
991 /*     estimate that is wrong by at most a factor of SQRT(2) */
992     spdiam = gu - gl;
993 /*     Compute splitting points */
994     dlarra_(n, &d__[1], &e[1], &e2[1], spltol, &spdiam, nsplit, &isplit[1], &
995             iinfo);
996 /*     Can force use of bisection instead of faster DQDS. */
997 /*     Option left in the code for future multisection work. */
998     forceb = FALSE_;
999 /*     Initialize USEDQD, DQDS should be used for ALLRNG unless someone */
1000 /*     explicitly wants bisection. */
1001     usedqd = irange == 1 && ! forceb;
1002     if (irange == 1 && ! forceb) {
1003 /*        Set interval [VL,VU] that contains all eigenvalues */
1004         *vl = gl;
1005         *vu = gu;
1006     } else {
1007 /*        We call DLARRD to find crude approximations to the eigenvalues */
1008 /*        in the desired range. In case IRANGE = INDRNG, we also obtain the */
1009 /*        interval (VL,VU] that contains all the wanted eigenvalues. */
1010 /*        An interval [LEFT,RIGHT] has converged if */
1011 /*        RIGHT-LEFT.LT.RTOL*MAX(ABS(LEFT),ABS(RIGHT)) */
1012 /*        DLARRD needs a WORK of size 4*N, IWORK of size 3*N */
1013         dlarrd_(range, "B", n, vl, vu, il, iu, &gers[1], &bsrtol, &d__[1], &e[
1014                 1], &e2[1], pivmin, nsplit, &isplit[1], &mm, &w[1], &werr[1], 
1015                 vl, vu, &iblock[1], &indexw[1], &work[1], &iwork[1], &iinfo);
1016         if (iinfo != 0) {
1017             *info = -1;
1018             return 0;
1019         }
1020 /*        Make sure that the entries M+1 to N in W, WERR, IBLOCK, INDEXW are 0 */
1021         i__1 = *n;
1022         for (i__ = mm + 1; i__ <= i__1; ++i__) {
1023             w[i__] = 0.;
1024             werr[i__] = 0.;
1025             iblock[i__] = 0;
1026             indexw[i__] = 0;
1027 /* L14: */
1028         }
1029     }
1030 /* ** */
1031 /*     Loop over unreduced blocks */
1032     ibegin = 1;
1033     wbegin = 1;
1034     i__1 = *nsplit;
1035     for (jblk = 1; jblk <= i__1; ++jblk) {
1036         iend = isplit[jblk];
1037         in = iend - ibegin + 1;
1038 /*        1 X 1 block */
1039         if (in == 1) {
1040             if (irange == 1 || irange == 3 && d__[ibegin] > *vl && d__[ibegin]
1041                      <= *vu || irange == 2 && iblock[wbegin] == jblk) {
1042                 ++(*m);
1043                 w[*m] = d__[ibegin];
1044                 werr[*m] = 0.;
1045 /*              The gap for a single block doesn't matter for the later */
1046 /*              algorithm and is assigned an arbitrary large value */
1047                 wgap[*m] = 0.;
1048                 iblock[*m] = jblk;
1049                 indexw[*m] = 1;
1050                 ++wbegin;
1051             }
1052 /*           E( IEND ) holds the shift for the initial RRR */
1053             e[iend] = 0.;
1054             ibegin = iend + 1;
1055             goto L170;
1056         }
1057
1058 /*        Blocks of size larger than 1x1 */
1059
1060 /*        E( IEND ) will hold the shift for the initial RRR, for now set it =0 */
1061         e[iend] = 0.;
1062
1063 /*        Find local outer bounds GL,GU for the block */
1064         gl = d__[ibegin];
1065         gu = d__[ibegin];
1066         i__2 = iend;
1067         for (i__ = ibegin; i__ <= i__2; ++i__) {
1068 /* Computing MIN */
1069             d__1 = gers[(i__ << 1) - 1];
1070             gl = f2cmin(d__1,gl);
1071 /* Computing MAX */
1072             d__1 = gers[i__ * 2];
1073             gu = f2cmax(d__1,gu);
1074 /* L15: */
1075         }
1076         spdiam = gu - gl;
1077         if (! (irange == 1 && ! forceb)) {
1078 /*           Count the number of eigenvalues in the current block. */
1079             mb = 0;
1080             i__2 = mm;
1081             for (i__ = wbegin; i__ <= i__2; ++i__) {
1082                 if (iblock[i__] == jblk) {
1083                     ++mb;
1084                 } else {
1085                     goto L21;
1086                 }
1087 /* L20: */
1088             }
1089 L21:
1090             if (mb == 0) {
1091 /*              No eigenvalue in the current block lies in the desired range */
1092 /*              E( IEND ) holds the shift for the initial RRR */
1093                 e[iend] = 0.;
1094                 ibegin = iend + 1;
1095                 goto L170;
1096             } else {
1097 /*              Decide whether dqds or bisection is more efficient */
1098                 usedqd = (doublereal) mb > in * .5 && ! forceb;
1099                 wend = wbegin + mb - 1;
1100 /*              Calculate gaps for the current block */
1101 /*              In later stages, when representations for individual */
1102 /*              eigenvalues are different, we use SIGMA = E( IEND ). */
1103                 sigma = 0.;
1104                 i__2 = wend - 1;
1105                 for (i__ = wbegin; i__ <= i__2; ++i__) {
1106 /* Computing MAX */
1107                     d__1 = 0., d__2 = w[i__ + 1] - werr[i__ + 1] - (w[i__] + 
1108                             werr[i__]);
1109                     wgap[i__] = f2cmax(d__1,d__2);
1110 /* L30: */
1111                 }
1112 /* Computing MAX */
1113                 d__1 = 0., d__2 = *vu - sigma - (w[wend] + werr[wend]);
1114                 wgap[wend] = f2cmax(d__1,d__2);
1115 /*              Find local index of the first and last desired evalue. */
1116                 indl = indexw[wbegin];
1117                 indu = indexw[wend];
1118             }
1119         }
1120         if (irange == 1 && ! forceb || usedqd) {
1121 /*           Case of DQDS */
1122 /*           Find approximations to the extremal eigenvalues of the block */
1123             dlarrk_(&in, &c__1, &gl, &gu, &d__[ibegin], &e2[ibegin], pivmin, &
1124                     rtl, &tmp, &tmp1, &iinfo);
1125             if (iinfo != 0) {
1126                 *info = -1;
1127                 return 0;
1128             }
1129 /* Computing MAX */
1130             d__2 = gl, d__3 = tmp - tmp1 - eps * 100. * (d__1 = tmp - tmp1, 
1131                     abs(d__1));
1132             isleft = f2cmax(d__2,d__3);
1133             dlarrk_(&in, &in, &gl, &gu, &d__[ibegin], &e2[ibegin], pivmin, &
1134                     rtl, &tmp, &tmp1, &iinfo);
1135             if (iinfo != 0) {
1136                 *info = -1;
1137                 return 0;
1138             }
1139 /* Computing MIN */
1140             d__2 = gu, d__3 = tmp + tmp1 + eps * 100. * (d__1 = tmp + tmp1, 
1141                     abs(d__1));
1142             isrght = f2cmin(d__2,d__3);
1143 /*           Improve the estimate of the spectral diameter */
1144             spdiam = isrght - isleft;
1145         } else {
1146 /*           Case of bisection */
1147 /*           Find approximations to the wanted extremal eigenvalues */
1148 /* Computing MAX */
1149             d__2 = gl, d__3 = w[wbegin] - werr[wbegin] - eps * 100. * (d__1 = 
1150                     w[wbegin] - werr[wbegin], abs(d__1));
1151             isleft = f2cmax(d__2,d__3);
1152 /* Computing MIN */
1153             d__2 = gu, d__3 = w[wend] + werr[wend] + eps * 100. * (d__1 = w[
1154                     wend] + werr[wend], abs(d__1));
1155             isrght = f2cmin(d__2,d__3);
1156         }
1157 /*        Decide whether the base representation for the current block */
1158 /*        L_JBLK D_JBLK L_JBLK^T = T_JBLK - sigma_JBLK I */
1159 /*        should be on the left or the right end of the current block. */
1160 /*        The strategy is to shift to the end which is "more populated" */
1161 /*        Furthermore, decide whether to use DQDS for the computation of */
1162 /*        the eigenvalue approximations at the end of DLARRE or bisection. */
1163 /*        dqds is chosen if all eigenvalues are desired or the number of */
1164 /*        eigenvalues to be computed is large compared to the blocksize. */
1165         if (irange == 1 && ! forceb) {
1166 /*           If all the eigenvalues have to be computed, we use dqd */
1167             usedqd = TRUE_;
1168 /*           INDL is the local index of the first eigenvalue to compute */
1169             indl = 1;
1170             indu = in;
1171 /*           MB =  number of eigenvalues to compute */
1172             mb = in;
1173             wend = wbegin + mb - 1;
1174 /*           Define 1/4 and 3/4 points of the spectrum */
1175             s1 = isleft + spdiam * .25;
1176             s2 = isrght - spdiam * .25;
1177         } else {
1178 /*           DLARRD has computed IBLOCK and INDEXW for each eigenvalue */
1179 /*           approximation. */
1180 /*           choose sigma */
1181             if (usedqd) {
1182                 s1 = isleft + spdiam * .25;
1183                 s2 = isrght - spdiam * .25;
1184             } else {
1185                 tmp = f2cmin(isrght,*vu) - f2cmax(isleft,*vl);
1186                 s1 = f2cmax(isleft,*vl) + tmp * .25;
1187                 s2 = f2cmin(isrght,*vu) - tmp * .25;
1188             }
1189         }
1190 /*        Compute the negcount at the 1/4 and 3/4 points */
1191         if (mb > 1) {
1192             dlarrc_("T", &in, &s1, &s2, &d__[ibegin], &e[ibegin], pivmin, &
1193                     cnt, &cnt1, &cnt2, &iinfo);
1194         }
1195         if (mb == 1) {
1196             sigma = gl;
1197             sgndef = 1.;
1198         } else if (cnt1 - indl >= indu - cnt2) {
1199             if (irange == 1 && ! forceb) {
1200                 sigma = f2cmax(isleft,gl);
1201             } else if (usedqd) {
1202 /*              use Gerschgorin bound as shift to get pos def matrix */
1203 /*              for dqds */
1204                 sigma = isleft;
1205             } else {
1206 /*              use approximation of the first desired eigenvalue of the */
1207 /*              block as shift */
1208                 sigma = f2cmax(isleft,*vl);
1209             }
1210             sgndef = 1.;
1211         } else {
1212             if (irange == 1 && ! forceb) {
1213                 sigma = f2cmin(isrght,gu);
1214             } else if (usedqd) {
1215 /*              use Gerschgorin bound as shift to get neg def matrix */
1216 /*              for dqds */
1217                 sigma = isrght;
1218             } else {
1219 /*              use approximation of the first desired eigenvalue of the */
1220 /*              block as shift */
1221                 sigma = f2cmin(isrght,*vu);
1222             }
1223             sgndef = -1.;
1224         }
1225 /*        An initial SIGMA has been chosen that will be used for computing */
1226 /*        T - SIGMA I = L D L^T */
1227 /*        Define the increment TAU of the shift in case the initial shift */
1228 /*        needs to be refined to obtain a factorization with not too much */
1229 /*        element growth. */
1230         if (usedqd) {
1231 /*           The initial SIGMA was to the outer end of the spectrum */
1232 /*           the matrix is definite and we need not retreat. */
1233             tau = spdiam * eps * *n + *pivmin * 2.;
1234 /* Computing MAX */
1235             d__1 = tau, d__2 = eps * 2. * abs(sigma);
1236             tau = f2cmax(d__1,d__2);
1237         } else {
1238             if (mb > 1) {
1239                 clwdth = w[wend] + werr[wend] - w[wbegin] - werr[wbegin];
1240                 avgap = (d__1 = clwdth / (doublereal) (wend - wbegin), abs(
1241                         d__1));
1242                 if (sgndef == 1.) {
1243 /* Computing MAX */
1244                     d__1 = wgap[wbegin];
1245                     tau = f2cmax(d__1,avgap) * .5;
1246 /* Computing MAX */
1247                     d__1 = tau, d__2 = werr[wbegin];
1248                     tau = f2cmax(d__1,d__2);
1249                 } else {
1250 /* Computing MAX */
1251                     d__1 = wgap[wend - 1];
1252                     tau = f2cmax(d__1,avgap) * .5;
1253 /* Computing MAX */
1254                     d__1 = tau, d__2 = werr[wend];
1255                     tau = f2cmax(d__1,d__2);
1256                 }
1257             } else {
1258                 tau = werr[wbegin];
1259             }
1260         }
1261
1262         for (idum = 1; idum <= 6; ++idum) {
1263 /*           Compute L D L^T factorization of tridiagonal matrix T - sigma I. */
1264 /*           Store D in WORK(1:IN), L in WORK(IN+1:2*IN), and reciprocals of */
1265 /*           pivots in WORK(2*IN+1:3*IN) */
1266             dpivot = d__[ibegin] - sigma;
1267             work[1] = dpivot;
1268             dmax__ = abs(work[1]);
1269             j = ibegin;
1270             i__2 = in - 1;
1271             for (i__ = 1; i__ <= i__2; ++i__) {
1272                 work[(in << 1) + i__] = 1. / work[i__];
1273                 tmp = e[j] * work[(in << 1) + i__];
1274                 work[in + i__] = tmp;
1275                 dpivot = d__[j + 1] - sigma - tmp * e[j];
1276                 work[i__ + 1] = dpivot;
1277 /* Computing MAX */
1278                 d__1 = dmax__, d__2 = abs(dpivot);
1279                 dmax__ = f2cmax(d__1,d__2);
1280                 ++j;
1281 /* L70: */
1282             }
1283 /*           check for element growth */
1284             if (dmax__ > spdiam * 64.) {
1285                 norep = TRUE_;
1286             } else {
1287                 norep = FALSE_;
1288             }
1289             if (usedqd && ! norep) {
1290 /*              Ensure the definiteness of the representation */
1291 /*              All entries of D (of L D L^T) must have the same sign */
1292                 i__2 = in;
1293                 for (i__ = 1; i__ <= i__2; ++i__) {
1294                     tmp = sgndef * work[i__];
1295                     if (tmp < 0.) {
1296                         norep = TRUE_;
1297                     }
1298 /* L71: */
1299                 }
1300             }
1301             if (norep) {
1302 /*              Note that in the case of IRANGE=ALLRNG, we use the Gerschgorin */
1303 /*              shift which makes the matrix definite. So we should end up */
1304 /*              here really only in the case of IRANGE = VALRNG or INDRNG. */
1305                 if (idum == 5) {
1306                     if (sgndef == 1.) {
1307 /*                    The fudged Gerschgorin shift should succeed */
1308                         sigma = gl - spdiam * 2. * eps * *n - *pivmin * 4.;
1309                     } else {
1310                         sigma = gu + spdiam * 2. * eps * *n + *pivmin * 4.;
1311                     }
1312                 } else {
1313                     sigma -= sgndef * tau;
1314                     tau *= 2.;
1315                 }
1316             } else {
1317 /*              an initial RRR is found */
1318                 goto L83;
1319             }
1320 /* L80: */
1321         }
1322 /*        if the program reaches this point, no base representation could be */
1323 /*        found in MAXTRY iterations. */
1324         *info = 2;
1325         return 0;
1326 L83:
1327 /*        At this point, we have found an initial base representation */
1328 /*        T - SIGMA I = L D L^T with not too much element growth. */
1329 /*        Store the shift. */
1330         e[iend] = sigma;
1331 /*        Store D and L. */
1332         dcopy_(&in, &work[1], &c__1, &d__[ibegin], &c__1);
1333         i__2 = in - 1;
1334         dcopy_(&i__2, &work[in + 1], &c__1, &e[ibegin], &c__1);
1335         if (mb > 1) {
1336
1337 /*           Perturb each entry of the base representation by a small */
1338 /*           (but random) relative amount to overcome difficulties with */
1339 /*           glued matrices. */
1340
1341             for (i__ = 1; i__ <= 4; ++i__) {
1342                 iseed[i__ - 1] = 1;
1343 /* L122: */
1344             }
1345             i__2 = (in << 1) - 1;
1346             dlarnv_(&c__2, iseed, &i__2, &work[1]);
1347             i__2 = in - 1;
1348             for (i__ = 1; i__ <= i__2; ++i__) {
1349                 d__[ibegin + i__ - 1] *= eps * 8. * work[i__] + 1.;
1350                 e[ibegin + i__ - 1] *= eps * 8. * work[in + i__] + 1.;
1351 /* L125: */
1352             }
1353             d__[iend] *= eps * 4. * work[in] + 1.;
1354
1355         }
1356
1357 /*        Don't update the Gerschgorin intervals because keeping track */
1358 /*        of the updates would be too much work in DLARRV. */
1359 /*        We update W instead and use it to locate the proper Gerschgorin */
1360 /*        intervals. */
1361 /*        Compute the required eigenvalues of L D L' by bisection or dqds */
1362         if (! usedqd) {
1363 /*           If DLARRD has been used, shift the eigenvalue approximations */
1364 /*           according to their representation. This is necessary for */
1365 /*           a uniform DLARRV since dqds computes eigenvalues of the */
1366 /*           shifted representation. In DLARRV, W will always hold the */
1367 /*           UNshifted eigenvalue approximation. */
1368             i__2 = wend;
1369             for (j = wbegin; j <= i__2; ++j) {
1370                 w[j] -= sigma;
1371                 werr[j] += (d__1 = w[j], abs(d__1)) * eps;
1372 /* L134: */
1373             }
1374 /*           call DLARRB to reduce eigenvalue error of the approximations */
1375 /*           from DLARRD */
1376             i__2 = iend - 1;
1377             for (i__ = ibegin; i__ <= i__2; ++i__) {
1378 /* Computing 2nd power */
1379                 d__1 = e[i__];
1380                 work[i__] = d__[i__] * (d__1 * d__1);
1381 /* L135: */
1382             }
1383 /*           use bisection to find EV from INDL to INDU */
1384             i__2 = indl - 1;
1385             dlarrb_(&in, &d__[ibegin], &work[ibegin], &indl, &indu, rtol1, 
1386                     rtol2, &i__2, &w[wbegin], &wgap[wbegin], &werr[wbegin], &
1387                     work[(*n << 1) + 1], &iwork[1], pivmin, &spdiam, &in, &
1388                     iinfo);
1389             if (iinfo != 0) {
1390                 *info = -4;
1391                 return 0;
1392             }
1393 /*           DLARRB computes all gaps correctly except for the last one */
1394 /*           Record distance to VU/GU */
1395 /* Computing MAX */
1396             d__1 = 0., d__2 = *vu - sigma - (w[wend] + werr[wend]);
1397             wgap[wend] = f2cmax(d__1,d__2);
1398             i__2 = indu;
1399             for (i__ = indl; i__ <= i__2; ++i__) {
1400                 ++(*m);
1401                 iblock[*m] = jblk;
1402                 indexw[*m] = i__;
1403 /* L138: */
1404             }
1405         } else {
1406 /*           Call dqds to get all eigs (and then possibly delete unwanted */
1407 /*           eigenvalues). */
1408 /*           Note that dqds finds the eigenvalues of the L D L^T representation */
1409 /*           of T to high relative accuracy. High relative accuracy */
1410 /*           might be lost when the shift of the RRR is subtracted to obtain */
1411 /*           the eigenvalues of T. However, T is not guaranteed to define its */
1412 /*           eigenvalues to high relative accuracy anyway. */
1413 /*           Set RTOL to the order of the tolerance used in DLASQ2 */
1414 /*           This is an ESTIMATED error, the worst case bound is 4*N*EPS */
1415 /*           which is usually too large and requires unnecessary work to be */
1416 /*           done by bisection when computing the eigenvectors */
1417             rtol = log((doublereal) in) * 4. * eps;
1418             j = ibegin;
1419             i__2 = in - 1;
1420             for (i__ = 1; i__ <= i__2; ++i__) {
1421                 work[(i__ << 1) - 1] = (d__1 = d__[j], abs(d__1));
1422                 work[i__ * 2] = e[j] * e[j] * work[(i__ << 1) - 1];
1423                 ++j;
1424 /* L140: */
1425             }
1426             work[(in << 1) - 1] = (d__1 = d__[iend], abs(d__1));
1427             work[in * 2] = 0.;
1428             dlasq2_(&in, &work[1], &iinfo);
1429             if (iinfo != 0) {
1430 /*              If IINFO = -5 then an index is part of a tight cluster */
1431 /*              and should be changed. The index is in IWORK(1) and the */
1432 /*              gap is in WORK(N+1) */
1433                 *info = -5;
1434                 return 0;
1435             } else {
1436 /*              Test that all eigenvalues are positive as expected */
1437                 i__2 = in;
1438                 for (i__ = 1; i__ <= i__2; ++i__) {
1439                     if (work[i__] < 0.) {
1440                         *info = -6;
1441                         return 0;
1442                     }
1443 /* L149: */
1444                 }
1445             }
1446             if (sgndef > 0.) {
1447                 i__2 = indu;
1448                 for (i__ = indl; i__ <= i__2; ++i__) {
1449                     ++(*m);
1450                     w[*m] = work[in - i__ + 1];
1451                     iblock[*m] = jblk;
1452                     indexw[*m] = i__;
1453 /* L150: */
1454                 }
1455             } else {
1456                 i__2 = indu;
1457                 for (i__ = indl; i__ <= i__2; ++i__) {
1458                     ++(*m);
1459                     w[*m] = -work[i__];
1460                     iblock[*m] = jblk;
1461                     indexw[*m] = i__;
1462 /* L160: */
1463                 }
1464             }
1465             i__2 = *m;
1466             for (i__ = *m - mb + 1; i__ <= i__2; ++i__) {
1467 /*              the value of RTOL below should be the tolerance in DLASQ2 */
1468                 werr[i__] = rtol * (d__1 = w[i__], abs(d__1));
1469 /* L165: */
1470             }
1471             i__2 = *m - 1;
1472             for (i__ = *m - mb + 1; i__ <= i__2; ++i__) {
1473 /*              compute the right gap between the intervals */
1474 /* Computing MAX */
1475                 d__1 = 0., d__2 = w[i__ + 1] - werr[i__ + 1] - (w[i__] + werr[
1476                         i__]);
1477                 wgap[i__] = f2cmax(d__1,d__2);
1478 /* L166: */
1479             }
1480 /* Computing MAX */
1481             d__1 = 0., d__2 = *vu - sigma - (w[*m] + werr[*m]);
1482             wgap[*m] = f2cmax(d__1,d__2);
1483         }
1484 /*        proceed with next block */
1485         ibegin = iend + 1;
1486         wbegin = wend + 1;
1487 L170:
1488         ;
1489     }
1490
1491     return 0;
1492
1493 /*     end of DLARRE */
1494
1495 } /* dlarre_ */
1496