C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / dlaqr5.c
1 #include <math.h>
2 #include <stdlib.h>
3 #include <string.h>
4 #include <stdio.h>
5 #include <complex.h>
6 #ifdef complex
7 #undef complex
8 #endif
9 #ifdef I
10 #undef I
11 #endif
12
13 #if defined(_WIN64)
14 typedef long long BLASLONG;
15 typedef unsigned long long BLASULONG;
16 #else
17 typedef long BLASLONG;
18 typedef unsigned long BLASULONG;
19 #endif
20
21 #ifdef LAPACK_ILP64
22 typedef BLASLONG blasint;
23 #if defined(_WIN64)
24 #define blasabs(x) llabs(x)
25 #else
26 #define blasabs(x) labs(x)
27 #endif
28 #else
29 typedef int blasint;
30 #define blasabs(x) abs(x)
31 #endif
32
33 typedef blasint integer;
34
35 typedef unsigned int uinteger;
36 typedef char *address;
37 typedef short int shortint;
38 typedef float real;
39 typedef double doublereal;
40 typedef struct { real r, i; } complex;
41 typedef struct { doublereal r, i; } doublecomplex;
42 #ifdef _MSC_VER
43 static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
44 static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
45 static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
46 static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
47 #else
48 static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
49 static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
50 static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
51 static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
52 #endif
53 #define pCf(z) (*_pCf(z))
54 #define pCd(z) (*_pCd(z))
55 typedef int logical;
56 typedef short int shortlogical;
57 typedef char logical1;
58 typedef char integer1;
59
60 #define TRUE_ (1)
61 #define FALSE_ (0)
62
63 /* Extern is for use with -E */
64 #ifndef Extern
65 #define Extern extern
66 #endif
67
68 /* I/O stuff */
69
70 typedef int flag;
71 typedef int ftnlen;
72 typedef int ftnint;
73
74 /*external read, write*/
75 typedef struct
76 {       flag cierr;
77         ftnint ciunit;
78         flag ciend;
79         char *cifmt;
80         ftnint cirec;
81 } cilist;
82
83 /*internal read, write*/
84 typedef struct
85 {       flag icierr;
86         char *iciunit;
87         flag iciend;
88         char *icifmt;
89         ftnint icirlen;
90         ftnint icirnum;
91 } icilist;
92
93 /*open*/
94 typedef struct
95 {       flag oerr;
96         ftnint ounit;
97         char *ofnm;
98         ftnlen ofnmlen;
99         char *osta;
100         char *oacc;
101         char *ofm;
102         ftnint orl;
103         char *oblnk;
104 } olist;
105
106 /*close*/
107 typedef struct
108 {       flag cerr;
109         ftnint cunit;
110         char *csta;
111 } cllist;
112
113 /*rewind, backspace, endfile*/
114 typedef struct
115 {       flag aerr;
116         ftnint aunit;
117 } alist;
118
119 /* inquire */
120 typedef struct
121 {       flag inerr;
122         ftnint inunit;
123         char *infile;
124         ftnlen infilen;
125         ftnint  *inex;  /*parameters in standard's order*/
126         ftnint  *inopen;
127         ftnint  *innum;
128         ftnint  *innamed;
129         char    *inname;
130         ftnlen  innamlen;
131         char    *inacc;
132         ftnlen  inacclen;
133         char    *inseq;
134         ftnlen  inseqlen;
135         char    *indir;
136         ftnlen  indirlen;
137         char    *infmt;
138         ftnlen  infmtlen;
139         char    *inform;
140         ftnint  informlen;
141         char    *inunf;
142         ftnlen  inunflen;
143         ftnint  *inrecl;
144         ftnint  *innrec;
145         char    *inblank;
146         ftnlen  inblanklen;
147 } inlist;
148
149 #define VOID void
150
151 union Multitype {       /* for multiple entry points */
152         integer1 g;
153         shortint h;
154         integer i;
155         /* longint j; */
156         real r;
157         doublereal d;
158         complex c;
159         doublecomplex z;
160         };
161
162 typedef union Multitype Multitype;
163
164 struct Vardesc {        /* for Namelist */
165         char *name;
166         char *addr;
167         ftnlen *dims;
168         int  type;
169         };
170 typedef struct Vardesc Vardesc;
171
172 struct Namelist {
173         char *name;
174         Vardesc **vars;
175         int nvars;
176         };
177 typedef struct Namelist Namelist;
178
179 #define abs(x) ((x) >= 0 ? (x) : -(x))
180 #define dabs(x) (fabs(x))
181 #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
182 #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
183 #define dmin(a,b) (f2cmin(a,b))
184 #define dmax(a,b) (f2cmax(a,b))
185 #define bit_test(a,b)   ((a) >> (b) & 1)
186 #define bit_clear(a,b)  ((a) & ~((uinteger)1 << (b)))
187 #define bit_set(a,b)    ((a) |  ((uinteger)1 << (b)))
188
189 #define abort_() { sig_die("Fortran abort routine called", 1); }
190 #define c_abs(z) (cabsf(Cf(z)))
191 #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
192 #ifdef _MSC_VER
193 #define c_div(c, a, b) {Cf(c)._Val[0] = (Cf(a)._Val[0]/Cf(b)._Val[0]); Cf(c)._Val[1]=(Cf(a)._Val[1]/Cf(b)._Val[1]);}
194 #define z_div(c, a, b) {Cd(c)._Val[0] = (Cd(a)._Val[0]/Cd(b)._Val[0]); Cd(c)._Val[1]=(Cd(a)._Val[1]/Cd(b)._Val[1]);}
195 #else
196 #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
197 #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
198 #endif
199 #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
200 #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
201 #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
202 //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
203 #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
204 #define d_abs(x) (fabs(*(x)))
205 #define d_acos(x) (acos(*(x)))
206 #define d_asin(x) (asin(*(x)))
207 #define d_atan(x) (atan(*(x)))
208 #define d_atn2(x, y) (atan2(*(x),*(y)))
209 #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
210 #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
211 #define d_cos(x) (cos(*(x)))
212 #define d_cosh(x) (cosh(*(x)))
213 #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
214 #define d_exp(x) (exp(*(x)))
215 #define d_imag(z) (cimag(Cd(z)))
216 #define r_imag(z) (cimagf(Cf(z)))
217 #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
218 #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
219 #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
220 #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
221 #define d_log(x) (log(*(x)))
222 #define d_mod(x, y) (fmod(*(x), *(y)))
223 #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
224 #define d_nint(x) u_nint(*(x))
225 #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
226 #define d_sign(a,b) u_sign(*(a),*(b))
227 #define r_sign(a,b) u_sign(*(a),*(b))
228 #define d_sin(x) (sin(*(x)))
229 #define d_sinh(x) (sinh(*(x)))
230 #define d_sqrt(x) (sqrt(*(x)))
231 #define d_tan(x) (tan(*(x)))
232 #define d_tanh(x) (tanh(*(x)))
233 #define i_abs(x) abs(*(x))
234 #define i_dnnt(x) ((integer)u_nint(*(x)))
235 #define i_len(s, n) (n)
236 #define i_nint(x) ((integer)u_nint(*(x)))
237 #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
238 #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
239 #define pow_si(B,E) spow_ui(*(B),*(E))
240 #define pow_ri(B,E) spow_ui(*(B),*(E))
241 #define pow_di(B,E) dpow_ui(*(B),*(E))
242 #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
243 #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
244 #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
245 #define s_cat(lpp, rpp, rnp, np, llp) {         ftnlen i, nc, ll; char *f__rp, *lp;     ll = (llp); lp = (lpp);         for(i=0; i < (int)*(np); ++i) {                 nc = ll;                if((rnp)[i] < nc) nc = (rnp)[i];                ll -= nc;               f__rp = (rpp)[i];               while(--nc >= 0) *lp++ = *(f__rp)++;         }  while(--ll >= 0) *lp++ = ' '; }
246 #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
247 #define s_copy(A,B,C,D) { int __i,__m; for (__i=0, __m=f2cmin((C),(D)); __i<__m && (B)[__i] != 0; ++__i) (A)[__i] = (B)[__i]; }
248 #define sig_die(s, kill) { exit(1); }
249 #define s_stop(s, n) {exit(0);}
250 static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
251 #define z_abs(z) (cabs(Cd(z)))
252 #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
253 #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
254 #define myexit_() break;
255 #define mycycle_() continue;
256 #define myceiling_(w) {ceil(w)}
257 #define myhuge_(w) {HUGE_VAL}
258 //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
259 #define mymaxloc_(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
260
261 /* procedure parameter types for -A and -C++ */
262
263 #define F2C_proc_par_types 1
264 #ifdef __cplusplus
265 typedef logical (*L_fp)(...);
266 #else
267 typedef logical (*L_fp)();
268 #endif
269
270 static float spow_ui(float x, integer n) {
271         float pow=1.0; unsigned long int u;
272         if(n != 0) {
273                 if(n < 0) n = -n, x = 1/x;
274                 for(u = n; ; ) {
275                         if(u & 01) pow *= x;
276                         if(u >>= 1) x *= x;
277                         else break;
278                 }
279         }
280         return pow;
281 }
282 static double dpow_ui(double x, integer n) {
283         double pow=1.0; unsigned long int u;
284         if(n != 0) {
285                 if(n < 0) n = -n, x = 1/x;
286                 for(u = n; ; ) {
287                         if(u & 01) pow *= x;
288                         if(u >>= 1) x *= x;
289                         else break;
290                 }
291         }
292         return pow;
293 }
294 #ifdef _MSC_VER
295 static _Fcomplex cpow_ui(complex x, integer n) {
296         complex pow={1.0,0.0}; unsigned long int u;
297                 if(n != 0) {
298                 if(n < 0) n = -n, x.r = 1/x.r, x.i=1/x.i;
299                 for(u = n; ; ) {
300                         if(u & 01) pow.r *= x.r, pow.i *= x.i;
301                         if(u >>= 1) x.r *= x.r, x.i *= x.i;
302                         else break;
303                 }
304         }
305         _Fcomplex p={pow.r, pow.i};
306         return p;
307 }
308 #else
309 static _Complex float cpow_ui(_Complex float x, integer n) {
310         _Complex float pow=1.0; unsigned long int u;
311         if(n != 0) {
312                 if(n < 0) n = -n, x = 1/x;
313                 for(u = n; ; ) {
314                         if(u & 01) pow *= x;
315                         if(u >>= 1) x *= x;
316                         else break;
317                 }
318         }
319         return pow;
320 }
321 #endif
322 #ifdef _MSC_VER
323 static _Dcomplex zpow_ui(_Dcomplex x, integer n) {
324         _Dcomplex pow={1.0,0.0}; unsigned long int u;
325         if(n != 0) {
326                 if(n < 0) n = -n, x._Val[0] = 1/x._Val[0], x._Val[1] =1/x._Val[1];
327                 for(u = n; ; ) {
328                         if(u & 01) pow._Val[0] *= x._Val[0], pow._Val[1] *= x._Val[1];
329                         if(u >>= 1) x._Val[0] *= x._Val[0], x._Val[1] *= x._Val[1];
330                         else break;
331                 }
332         }
333         _Dcomplex p = {pow._Val[0], pow._Val[1]};
334         return p;
335 }
336 #else
337 static _Complex double zpow_ui(_Complex double x, integer n) {
338         _Complex double pow=1.0; unsigned long int u;
339         if(n != 0) {
340                 if(n < 0) n = -n, x = 1/x;
341                 for(u = n; ; ) {
342                         if(u & 01) pow *= x;
343                         if(u >>= 1) x *= x;
344                         else break;
345                 }
346         }
347         return pow;
348 }
349 #endif
350 static integer pow_ii(integer x, integer n) {
351         integer pow; unsigned long int u;
352         if (n <= 0) {
353                 if (n == 0 || x == 1) pow = 1;
354                 else if (x != -1) pow = x == 0 ? 1/x : 0;
355                 else n = -n;
356         }
357         if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
358                 u = n;
359                 for(pow = 1; ; ) {
360                         if(u & 01) pow *= x;
361                         if(u >>= 1) x *= x;
362                         else break;
363                 }
364         }
365         return pow;
366 }
367 static integer dmaxloc_(double *w, integer s, integer e, integer *n)
368 {
369         double m; integer i, mi;
370         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
371                 if (w[i-1]>m) mi=i ,m=w[i-1];
372         return mi-s+1;
373 }
374 static integer smaxloc_(float *w, integer s, integer e, integer *n)
375 {
376         float m; integer i, mi;
377         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
378                 if (w[i-1]>m) mi=i ,m=w[i-1];
379         return mi-s+1;
380 }
381 static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
382         integer n = *n_, incx = *incx_, incy = *incy_, i;
383 #ifdef _MSC_VER
384         _Fcomplex zdotc = {0.0, 0.0};
385         if (incx == 1 && incy == 1) {
386                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
387                         zdotc._Val[0] += conjf(Cf(&x[i]))._Val[0] * Cf(&y[i])._Val[0];
388                         zdotc._Val[1] += conjf(Cf(&x[i]))._Val[1] * Cf(&y[i])._Val[1];
389                 }
390         } else {
391                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
392                         zdotc._Val[0] += conjf(Cf(&x[i*incx]))._Val[0] * Cf(&y[i*incy])._Val[0];
393                         zdotc._Val[1] += conjf(Cf(&x[i*incx]))._Val[1] * Cf(&y[i*incy])._Val[1];
394                 }
395         }
396         pCf(z) = zdotc;
397 }
398 #else
399         _Complex float zdotc = 0.0;
400         if (incx == 1 && incy == 1) {
401                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
402                         zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
403                 }
404         } else {
405                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
406                         zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
407                 }
408         }
409         pCf(z) = zdotc;
410 }
411 #endif
412 static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
413         integer n = *n_, incx = *incx_, incy = *incy_, i;
414 #ifdef _MSC_VER
415         _Dcomplex zdotc = {0.0, 0.0};
416         if (incx == 1 && incy == 1) {
417                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
418                         zdotc._Val[0] += conj(Cd(&x[i]))._Val[0] * Cd(&y[i])._Val[0];
419                         zdotc._Val[1] += conj(Cd(&x[i]))._Val[1] * Cd(&y[i])._Val[1];
420                 }
421         } else {
422                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
423                         zdotc._Val[0] += conj(Cd(&x[i*incx]))._Val[0] * Cd(&y[i*incy])._Val[0];
424                         zdotc._Val[1] += conj(Cd(&x[i*incx]))._Val[1] * Cd(&y[i*incy])._Val[1];
425                 }
426         }
427         pCd(z) = zdotc;
428 }
429 #else
430         _Complex double zdotc = 0.0;
431         if (incx == 1 && incy == 1) {
432                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
433                         zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
434                 }
435         } else {
436                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
437                         zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
438                 }
439         }
440         pCd(z) = zdotc;
441 }
442 #endif  
443 static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
444         integer n = *n_, incx = *incx_, incy = *incy_, i;
445 #ifdef _MSC_VER
446         _Fcomplex zdotc = {0.0, 0.0};
447         if (incx == 1 && incy == 1) {
448                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
449                         zdotc._Val[0] += Cf(&x[i])._Val[0] * Cf(&y[i])._Val[0];
450                         zdotc._Val[1] += Cf(&x[i])._Val[1] * Cf(&y[i])._Val[1];
451                 }
452         } else {
453                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
454                         zdotc._Val[0] += Cf(&x[i*incx])._Val[0] * Cf(&y[i*incy])._Val[0];
455                         zdotc._Val[1] += Cf(&x[i*incx])._Val[1] * Cf(&y[i*incy])._Val[1];
456                 }
457         }
458         pCf(z) = zdotc;
459 }
460 #else
461         _Complex float zdotc = 0.0;
462         if (incx == 1 && incy == 1) {
463                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
464                         zdotc += Cf(&x[i]) * Cf(&y[i]);
465                 }
466         } else {
467                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
468                         zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
469                 }
470         }
471         pCf(z) = zdotc;
472 }
473 #endif
474 static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
475         integer n = *n_, incx = *incx_, incy = *incy_, i;
476 #ifdef _MSC_VER
477         _Dcomplex zdotc = {0.0, 0.0};
478         if (incx == 1 && incy == 1) {
479                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
480                         zdotc._Val[0] += Cd(&x[i])._Val[0] * Cd(&y[i])._Val[0];
481                         zdotc._Val[1] += Cd(&x[i])._Val[1] * Cd(&y[i])._Val[1];
482                 }
483         } else {
484                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
485                         zdotc._Val[0] += Cd(&x[i*incx])._Val[0] * Cd(&y[i*incy])._Val[0];
486                         zdotc._Val[1] += Cd(&x[i*incx])._Val[1] * Cd(&y[i*incy])._Val[1];
487                 }
488         }
489         pCd(z) = zdotc;
490 }
491 #else
492         _Complex double zdotc = 0.0;
493         if (incx == 1 && incy == 1) {
494                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
495                         zdotc += Cd(&x[i]) * Cd(&y[i]);
496                 }
497         } else {
498                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
499                         zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
500                 }
501         }
502         pCd(z) = zdotc;
503 }
504 #endif
505 /*  -- translated by f2c (version 20000121).
506    You must link the resulting object file with the libraries:
507         -lf2c -lm   (in that order)
508 */
509
510
511
512
513
514 /* Table of constant values */
515
516 static doublereal c_b7 = 0.;
517 static doublereal c_b8 = 1.;
518 static integer c__2 = 2;
519 static integer c__1 = 1;
520 static integer c__3 = 3;
521
522 /* > \brief \b DLAQR5 performs a single small-bulge multi-shift QR sweep. */
523
524 /*  =========== DOCUMENTATION =========== */
525
526 /* Online html documentation available at */
527 /*            http://www.netlib.org/lapack/explore-html/ */
528
529 /* > \htmlonly */
530 /* > Download DLAQR5 + dependencies */
531 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlaqr5.
532 f"> */
533 /* > [TGZ]</a> */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlaqr5.
535 f"> */
536 /* > [ZIP]</a> */
537 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlaqr5.
538 f"> */
539 /* > [TXT]</a> */
540 /* > \endhtmlonly */
541
542 /*  Definition: */
543 /*  =========== */
544
545 /*       SUBROUTINE DLAQR5( WANTT, WANTZ, KACC22, N, KTOP, KBOT, NSHFTS, */
546 /*                          SR, SI, H, LDH, ILOZ, IHIZ, Z, LDZ, V, LDV, U, */
547 /*                          LDU, NV, WV, LDWV, NH, WH, LDWH ) */
548
549 /*       INTEGER            IHIZ, ILOZ, KACC22, KBOT, KTOP, LDH, LDU, LDV, */
550 /*      $                   LDWH, LDWV, LDZ, N, NH, NSHFTS, NV */
551 /*       LOGICAL            WANTT, WANTZ */
552 /*       DOUBLE PRECISION   H( LDH, * ), SI( * ), SR( * ), U( LDU, * ), */
553 /*      $                   V( LDV, * ), WH( LDWH, * ), WV( LDWV, * ), */
554 /*      $                   Z( LDZ, * ) */
555
556
557 /* > \par Purpose: */
558 /*  ============= */
559 /* > */
560 /* > \verbatim */
561 /* > */
562 /* >    DLAQR5, called by DLAQR0, performs a */
563 /* >    single small-bulge multi-shift QR sweep. */
564 /* > \endverbatim */
565
566 /*  Arguments: */
567 /*  ========== */
568
569 /* > \param[in] WANTT */
570 /* > \verbatim */
571 /* >          WANTT is LOGICAL */
572 /* >             WANTT = .true. if the quasi-triangular Schur factor */
573 /* >             is being computed.  WANTT is set to .false. otherwise. */
574 /* > \endverbatim */
575 /* > */
576 /* > \param[in] WANTZ */
577 /* > \verbatim */
578 /* >          WANTZ is LOGICAL */
579 /* >             WANTZ = .true. if the orthogonal Schur factor is being */
580 /* >             computed.  WANTZ is set to .false. otherwise. */
581 /* > \endverbatim */
582 /* > */
583 /* > \param[in] KACC22 */
584 /* > \verbatim */
585 /* >          KACC22 is INTEGER with value 0, 1, or 2. */
586 /* >             Specifies the computation mode of far-from-diagonal */
587 /* >             orthogonal updates. */
588 /* >        = 0: DLAQR5 does not accumulate reflections and does not */
589 /* >             use matrix-matrix multiply to update far-from-diagonal */
590 /* >             matrix entries. */
591 /* >        = 1: DLAQR5 accumulates reflections and uses matrix-matrix */
592 /* >             multiply to update the far-from-diagonal matrix entries. */
593 /* >        = 2: Same as KACC22 = 1. This option used to enable exploiting */
594 /* >             the 2-by-2 structure during matrix multiplications, but */
595 /* >             this is no longer supported. */
596 /* > \endverbatim */
597 /* > */
598 /* > \param[in] N */
599 /* > \verbatim */
600 /* >          N is INTEGER */
601 /* >             N is the order of the Hessenberg matrix H upon which this */
602 /* >             subroutine operates. */
603 /* > \endverbatim */
604 /* > */
605 /* > \param[in] KTOP */
606 /* > \verbatim */
607 /* >          KTOP is INTEGER */
608 /* > \endverbatim */
609 /* > */
610 /* > \param[in] KBOT */
611 /* > \verbatim */
612 /* >          KBOT is INTEGER */
613 /* >             These are the first and last rows and columns of an */
614 /* >             isolated diagonal block upon which the QR sweep is to be */
615 /* >             applied. It is assumed without a check that */
616 /* >                       either KTOP = 1  or   H(KTOP,KTOP-1) = 0 */
617 /* >             and */
618 /* >                       either KBOT = N  or   H(KBOT+1,KBOT) = 0. */
619 /* > \endverbatim */
620 /* > */
621 /* > \param[in] NSHFTS */
622 /* > \verbatim */
623 /* >          NSHFTS is INTEGER */
624 /* >             NSHFTS gives the number of simultaneous shifts.  NSHFTS */
625 /* >             must be positive and even. */
626 /* > \endverbatim */
627 /* > */
628 /* > \param[in,out] SR */
629 /* > \verbatim */
630 /* >          SR is DOUBLE PRECISION array, dimension (NSHFTS) */
631 /* > \endverbatim */
632 /* > */
633 /* > \param[in,out] SI */
634 /* > \verbatim */
635 /* >          SI is DOUBLE PRECISION array, dimension (NSHFTS) */
636 /* >             SR contains the real parts and SI contains the imaginary */
637 /* >             parts of the NSHFTS shifts of origin that define the */
638 /* >             multi-shift QR sweep.  On output SR and SI may be */
639 /* >             reordered. */
640 /* > \endverbatim */
641 /* > */
642 /* > \param[in,out] H */
643 /* > \verbatim */
644 /* >          H is DOUBLE PRECISION array, dimension (LDH,N) */
645 /* >             On input H contains a Hessenberg matrix.  On output a */
646 /* >             multi-shift QR sweep with shifts SR(J)+i*SI(J) is applied */
647 /* >             to the isolated diagonal block in rows and columns KTOP */
648 /* >             through KBOT. */
649 /* > \endverbatim */
650 /* > */
651 /* > \param[in] LDH */
652 /* > \verbatim */
653 /* >          LDH is INTEGER */
654 /* >             LDH is the leading dimension of H just as declared in the */
655 /* >             calling procedure.  LDH >= MAX(1,N). */
656 /* > \endverbatim */
657 /* > */
658 /* > \param[in] ILOZ */
659 /* > \verbatim */
660 /* >          ILOZ is INTEGER */
661 /* > \endverbatim */
662 /* > */
663 /* > \param[in] IHIZ */
664 /* > \verbatim */
665 /* >          IHIZ is INTEGER */
666 /* >             Specify the rows of Z to which transformations must be */
667 /* >             applied if WANTZ is .TRUE.. 1 <= ILOZ <= IHIZ <= N */
668 /* > \endverbatim */
669 /* > */
670 /* > \param[in,out] Z */
671 /* > \verbatim */
672 /* >          Z is DOUBLE PRECISION array, dimension (LDZ,IHIZ) */
673 /* >             If WANTZ = .TRUE., then the QR Sweep orthogonal */
674 /* >             similarity transformation is accumulated into */
675 /* >             Z(ILOZ:IHIZ,ILOZ:IHIZ) from the right. */
676 /* >             If WANTZ = .FALSE., then Z is unreferenced. */
677 /* > \endverbatim */
678 /* > */
679 /* > \param[in] LDZ */
680 /* > \verbatim */
681 /* >          LDZ is INTEGER */
682 /* >             LDA is the leading dimension of Z just as declared in */
683 /* >             the calling procedure. LDZ >= N. */
684 /* > \endverbatim */
685 /* > */
686 /* > \param[out] V */
687 /* > \verbatim */
688 /* >          V is DOUBLE PRECISION array, dimension (LDV,NSHFTS/2) */
689 /* > \endverbatim */
690 /* > */
691 /* > \param[in] LDV */
692 /* > \verbatim */
693 /* >          LDV is INTEGER */
694 /* >             LDV is the leading dimension of V as declared in the */
695 /* >             calling procedure.  LDV >= 3. */
696 /* > \endverbatim */
697 /* > */
698 /* > \param[out] U */
699 /* > \verbatim */
700 /* >          U is DOUBLE PRECISION array, dimension (LDU,2*NSHFTS) */
701 /* > \endverbatim */
702 /* > */
703 /* > \param[in] LDU */
704 /* > \verbatim */
705 /* >          LDU is INTEGER */
706 /* >             LDU is the leading dimension of U just as declared in the */
707 /* >             in the calling subroutine.  LDU >= 2*NSHFTS. */
708 /* > \endverbatim */
709 /* > */
710 /* > \param[in] NV */
711 /* > \verbatim */
712 /* >          NV is INTEGER */
713 /* >             NV is the number of rows in WV agailable for workspace. */
714 /* >             NV >= 1. */
715 /* > \endverbatim */
716 /* > */
717 /* > \param[out] WV */
718 /* > \verbatim */
719 /* >          WV is DOUBLE PRECISION array, dimension (LDWV,2*NSHFTS) */
720 /* > \endverbatim */
721 /* > */
722 /* > \param[in] LDWV */
723 /* > \verbatim */
724 /* >          LDWV is INTEGER */
725 /* >             LDWV is the leading dimension of WV as declared in the */
726 /* >             in the calling subroutine.  LDWV >= NV. */
727 /* > \endverbatim */
728
729 /* > \param[in] NH */
730 /* > \verbatim */
731 /* >          NH is INTEGER */
732 /* >             NH is the number of columns in array WH available for */
733 /* >             workspace. NH >= 1. */
734 /* > \endverbatim */
735 /* > */
736 /* > \param[out] WH */
737 /* > \verbatim */
738 /* >          WH is DOUBLE PRECISION array, dimension (LDWH,NH) */
739 /* > \endverbatim */
740 /* > */
741 /* > \param[in] LDWH */
742 /* > \verbatim */
743 /* >          LDWH is INTEGER */
744 /* >             Leading dimension of WH just as declared in the */
745 /* >             calling procedure.  LDWH >= 2*NSHFTS. */
746 /* > \endverbatim */
747 /* > */
748 /*  Authors: */
749 /*  ======== */
750
751 /* > \author Univ. of Tennessee */
752 /* > \author Univ. of California Berkeley */
753 /* > \author Univ. of Colorado Denver */
754 /* > \author NAG Ltd. */
755
756 /* > \date January 2021 */
757
758 /* > \ingroup doubleOTHERauxiliary */
759
760 /* > \par Contributors: */
761 /*  ================== */
762 /* > */
763 /* >       Karen Braman and Ralph Byers, Department of Mathematics, */
764 /* >       University of Kansas, USA */
765 /* > */
766 /* >       Lars Karlsson, Daniel Kressner, and Bruno Lang */
767 /* > */
768 /* >       Thijs Steel, Department of Computer science, */
769 /* >       KU Leuven, Belgium */
770
771 /* > \par References: */
772 /*  ================ */
773 /* > */
774 /* >       K. Braman, R. Byers and R. Mathias, The Multi-Shift QR */
775 /* >       Algorithm Part I: Maintaining Well Focused Shifts, and Level 3 */
776 /* >       Performance, SIAM Journal of Matrix Analysis, volume 23, pages */
777 /* >       929--947, 2002. */
778 /* > */
779 /* >       Lars Karlsson, Daniel Kressner, and Bruno Lang, Optimally packed */
780 /* >       chains of bulges in multishift QR algorithms. */
781 /* >       ACM Trans. Math. Softw. 40, 2, Article 12 (February 2014). */
782 /* > */
783 /*  ===================================================================== */
784 /* Subroutine */ int dlaqr5_(logical *wantt, logical *wantz, integer *kacc22, 
785         integer *n, integer *ktop, integer *kbot, integer *nshfts, doublereal 
786         *sr, doublereal *si, doublereal *h__, integer *ldh, integer *iloz, 
787         integer *ihiz, doublereal *z__, integer *ldz, doublereal *v, integer *
788         ldv, doublereal *u, integer *ldu, integer *nv, doublereal *wv, 
789         integer *ldwv, integer *nh, doublereal *wh, integer *ldwh)
790 {
791     /* System generated locals */
792     integer h_dim1, h_offset, u_dim1, u_offset, v_dim1, v_offset, wh_dim1, 
793             wh_offset, wv_dim1, wv_offset, z_dim1, z_offset, i__1, i__2, i__3,
794              i__4, i__5, i__6, i__7;
795     doublereal d__1, d__2, d__3, d__4, d__5;
796
797     /* Local variables */
798     doublereal beta;
799     logical bmp22;
800     integer jcol, jlen, jbot, mbot;
801     doublereal swap;
802     integer jtop, jrow, mtop, i__, j, k, m;
803     doublereal alpha;
804     logical accum;
805     extern /* Subroutine */ int dgemm_(char *, char *, integer *, integer *, 
806             integer *, doublereal *, doublereal *, integer *, doublereal *, 
807             integer *, doublereal *, doublereal *, integer *);
808     integer ndcol, incol, krcol, nbmps, i2, k1, i4;
809     extern /* Subroutine */ int dlaqr1_(integer *, doublereal *, integer *, 
810             doublereal *, doublereal *, doublereal *, doublereal *, 
811             doublereal *), dlabad_(doublereal *, doublereal *);
812     doublereal h11, h12, h21, h22;
813     integer m22;
814     extern doublereal dlamch_(char *);
815     extern /* Subroutine */ int dlarfg_(integer *, doublereal *, doublereal *,
816              integer *, doublereal *);
817     integer ns, nu;
818     doublereal vt[3];
819     extern /* Subroutine */ int dlacpy_(char *, integer *, integer *, 
820             doublereal *, integer *, doublereal *, integer *);
821     doublereal safmin, safmax;
822     extern /* Subroutine */ int dlaset_(char *, integer *, integer *, 
823             doublereal *, doublereal *, doublereal *, integer *);
824     doublereal refsum, smlnum, scl;
825     integer kdu, kms;
826     doublereal ulp;
827     doublereal tst1, tst2;
828
829
830 /*  -- LAPACK auxiliary routine (version 3.7.1) -- */
831 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
832 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
833 /*     June 2016 */
834
835
836 /*  ================================================================ */
837
838
839 /*     ==== If there are no shifts, then there is nothing to do. ==== */
840
841     /* Parameter adjustments */
842     --sr;
843     --si;
844     h_dim1 = *ldh;
845     h_offset = 1 + h_dim1 * 1;
846     h__ -= h_offset;
847     z_dim1 = *ldz;
848     z_offset = 1 + z_dim1 * 1;
849     z__ -= z_offset;
850     v_dim1 = *ldv;
851     v_offset = 1 + v_dim1 * 1;
852     v -= v_offset;
853     u_dim1 = *ldu;
854     u_offset = 1 + u_dim1 * 1;
855     u -= u_offset;
856     wv_dim1 = *ldwv;
857     wv_offset = 1 + wv_dim1 * 1;
858     wv -= wv_offset;
859     wh_dim1 = *ldwh;
860     wh_offset = 1 + wh_dim1 * 1;
861     wh -= wh_offset;
862
863     /* Function Body */
864     if (*nshfts < 2) {
865         return 0;
866     }
867
868 /*     ==== If the active block is empty or 1-by-1, then there */
869 /*     .    is nothing to do. ==== */
870
871     if (*ktop >= *kbot) {
872         return 0;
873     }
874
875 /*     ==== Shuffle shifts into pairs of real shifts and pairs */
876 /*     .    of complex conjugate shifts assuming complex */
877 /*     .    conjugate shifts are already adjacent to one */
878 /*     .    another. ==== */
879
880     i__1 = *nshfts - 2;
881     for (i__ = 1; i__ <= i__1; i__ += 2) {
882         if (si[i__] != -si[i__ + 1]) {
883
884             swap = sr[i__];
885             sr[i__] = sr[i__ + 1];
886             sr[i__ + 1] = sr[i__ + 2];
887             sr[i__ + 2] = swap;
888
889             swap = si[i__];
890             si[i__] = si[i__ + 1];
891             si[i__ + 1] = si[i__ + 2];
892             si[i__ + 2] = swap;
893         }
894 /* L10: */
895     }
896
897 /*     ==== NSHFTS is supposed to be even, but if it is odd, */
898 /*     .    then simply reduce it by one.  The shuffle above */
899 /*     .    ensures that the dropped shift is real and that */
900 /*     .    the remaining shifts are paired. ==== */
901
902     ns = *nshfts - *nshfts % 2;
903
904 /*     ==== Machine constants for deflation ==== */
905
906     safmin = dlamch_("SAFE MINIMUM");
907     safmax = 1. / safmin;
908     dlabad_(&safmin, &safmax);
909     ulp = dlamch_("PRECISION");
910     smlnum = safmin * ((doublereal) (*n) / ulp);
911
912 /*     ==== Use accumulated reflections to update far-from-diagonal */
913 /*     .    entries ? ==== */
914
915     accum = *kacc22 == 1 || *kacc22 == 2;
916
917 /*     ==== clear trash ==== */
918
919     if (*ktop + 2 <= *kbot) {
920         h__[*ktop + 2 + *ktop * h_dim1] = 0.;
921     }
922
923 /*     ==== NBMPS = number of 2-shift bulges in the chain ==== */
924
925     nbmps = ns / 2;
926
927 /*     ==== KDU = width of slab ==== */
928
929     kdu = nbmps << 2;
930
931 /*     ==== Create and chase chains of NBMPS bulges ==== */
932
933     i__1 = *kbot - 2;
934     i__2 = nbmps << 1;
935     for (incol = *ktop - (nbmps << 1) + 1; i__2 < 0 ? incol >= i__1 : incol <=
936              i__1; incol += i__2) {
937
938 /*        JTOP = Index from which updates from the right start. */
939
940         if (accum) {
941             jtop = f2cmax(*ktop,incol);
942         } else if (*wantt) {
943             jtop = 1;
944         } else {
945             jtop = *ktop;
946         }
947
948         ndcol = incol + kdu;
949         if (accum) {
950             dlaset_("ALL", &kdu, &kdu, &c_b7, &c_b8, &u[u_offset], ldu);
951         }
952
953 /*        ==== Near-the-diagonal bulge chase.  The following loop */
954 /*        .    performs the near-the-diagonal part of a small bulge */
955 /*        .    multi-shift QR sweep.  Each 4*NBMPS column diagonal */
956 /*        .    chunk extends from column INCOL to column NDCOL */
957 /*        .    (including both column INCOL and column NDCOL). The */
958 /*        .    following loop chases a 2*NBMPS+1 column long chain of */
959 /*        .    NBMPS bulges 2*NBMPS columns to the right.  (INCOL */
960 /*        .    may be less than KTOP and and NDCOL may be greater than */
961 /*        .    KBOT indicating phantom columns from which to chase */
962 /*        .    bulges before they are actually introduced or to which */
963 /*        .    to chase bulges beyond column KBOT.)  ==== */
964
965 /* Computing MIN */
966         i__4 = incol + (nbmps << 1) - 1, i__5 = *kbot - 2;
967         i__3 = f2cmin(i__4,i__5);
968         for (krcol = incol; krcol <= i__3; ++krcol) {
969
970 /*           ==== Bulges number MTOP to MBOT are active double implicit */
971 /*           .    shift bulges.  There may or may not also be small */
972 /*           .    2-by-2 bulge, if there is room.  The inactive bulges */
973 /*           .    (if any) must wait until the active bulges have moved */
974 /*           .    down the diagonal to make room.  The phantom matrix */
975 /*           .    paradigm described above helps keep track.  ==== */
976
977 /* Computing MAX */
978             i__4 = 1, i__5 = (*ktop - krcol) / 2 + 1;
979             mtop = f2cmax(i__4,i__5);
980 /* Computing MIN */
981             i__4 = nbmps, i__5 = (*kbot - krcol - 1) / 2;
982             mbot = f2cmin(i__4,i__5);
983             m22 = mbot + 1;
984             bmp22 = mbot < nbmps && krcol + (m22 - 1 << 1) == *kbot - 2;
985
986 /*           ==== Generate reflections to chase the chain right */
987 /*           .    one column.  (The minimum value of K is KTOP-1.) ==== */
988
989             if (bmp22) {
990
991 /*              ==== Special case: 2-by-2 reflection at bottom treated */
992 /*              .    separately ==== */
993
994                 k = krcol + (m22 - 1 << 1);
995                 if (k == *ktop - 1) {
996                     dlaqr1_(&c__2, &h__[k + 1 + (k + 1) * h_dim1], ldh, &sr[(
997                             m22 << 1) - 1], &si[(m22 << 1) - 1], &sr[m22 * 2],
998                              &si[m22 * 2], &v[m22 * v_dim1 + 1]);
999                     beta = v[m22 * v_dim1 + 1];
1000                     dlarfg_(&c__2, &beta, &v[m22 * v_dim1 + 2], &c__1, &v[m22 
1001                             * v_dim1 + 1]);
1002                 } else {
1003                     beta = h__[k + 1 + k * h_dim1];
1004                     v[m22 * v_dim1 + 2] = h__[k + 2 + k * h_dim1];
1005                     dlarfg_(&c__2, &beta, &v[m22 * v_dim1 + 2], &c__1, &v[m22 
1006                             * v_dim1 + 1]);
1007                     h__[k + 1 + k * h_dim1] = beta;
1008                     h__[k + 2 + k * h_dim1] = 0.;
1009                 }
1010
1011 /*              ==== Perform update from right within */
1012 /*              .    computational window. ==== */
1013
1014 /* Computing MIN */
1015                 i__5 = *kbot, i__6 = k + 3;
1016                 i__4 = f2cmin(i__5,i__6);
1017                 for (j = jtop; j <= i__4; ++j) {
1018                     refsum = v[m22 * v_dim1 + 1] * (h__[j + (k + 1) * h_dim1] 
1019                             + v[m22 * v_dim1 + 2] * h__[j + (k + 2) * h_dim1])
1020                             ;
1021                     h__[j + (k + 1) * h_dim1] -= refsum;
1022                     h__[j + (k + 2) * h_dim1] -= refsum * v[m22 * v_dim1 + 2];
1023 /* L30: */
1024                 }
1025
1026 /*              ==== Perform update from left within */
1027 /*              .    computational window. ==== */
1028
1029                 if (accum) {
1030                     jbot = f2cmin(ndcol,*kbot);
1031                 } else if (*wantt) {
1032                     jbot = *n;
1033                 } else {
1034                     jbot = *kbot;
1035                 }
1036                 i__4 = jbot;
1037                 for (j = k + 1; j <= i__4; ++j) {
1038                     refsum = v[m22 * v_dim1 + 1] * (h__[k + 1 + j * h_dim1] + 
1039                             v[m22 * v_dim1 + 2] * h__[k + 2 + j * h_dim1]);
1040                     h__[k + 1 + j * h_dim1] -= refsum;
1041                     h__[k + 2 + j * h_dim1] -= refsum * v[m22 * v_dim1 + 2];
1042 /* L40: */
1043                 }
1044
1045 /*              ==== The following convergence test requires that */
1046 /*              .    the tradition small-compared-to-nearby-diagonals */
1047 /*              .    criterion and the Ahues & Tisseur (LAWN 122, 1997) */
1048 /*              .    criteria both be satisfied.  The latter improves */
1049 /*              .    accuracy in some examples. Falling back on an */
1050 /*              .    alternate convergence criterion when TST1 or TST2 */
1051 /*              .    is zero (as done here) is traditional but probably */
1052 /*              .    unnecessary. ==== */
1053
1054                 if (k >= *ktop) {
1055                     if (h__[k + 1 + k * h_dim1] != 0.) {
1056                         tst1 = (d__1 = h__[k + k * h_dim1], abs(d__1)) + (
1057                                 d__2 = h__[k + 1 + (k + 1) * h_dim1], abs(
1058                                 d__2));
1059                         if (tst1 == 0.) {
1060                             if (k >= *ktop + 1) {
1061                                 tst1 += (d__1 = h__[k + (k - 1) * h_dim1], 
1062                                         abs(d__1));
1063                             }
1064                             if (k >= *ktop + 2) {
1065                                 tst1 += (d__1 = h__[k + (k - 2) * h_dim1], 
1066                                         abs(d__1));
1067                             }
1068                             if (k >= *ktop + 3) {
1069                                 tst1 += (d__1 = h__[k + (k - 3) * h_dim1], 
1070                                         abs(d__1));
1071                             }
1072                             if (k <= *kbot - 2) {
1073                                 tst1 += (d__1 = h__[k + 2 + (k + 1) * h_dim1],
1074                                          abs(d__1));
1075                             }
1076                             if (k <= *kbot - 3) {
1077                                 tst1 += (d__1 = h__[k + 3 + (k + 1) * h_dim1],
1078                                          abs(d__1));
1079                             }
1080                             if (k <= *kbot - 4) {
1081                                 tst1 += (d__1 = h__[k + 4 + (k + 1) * h_dim1],
1082                                          abs(d__1));
1083                             }
1084                         }
1085 /* Computing MAX */
1086                         d__2 = smlnum, d__3 = ulp * tst1;
1087                         if ((d__1 = h__[k + 1 + k * h_dim1], abs(d__1)) <= 
1088                                 f2cmax(d__2,d__3)) {
1089 /* Computing MAX */
1090                             d__3 = (d__1 = h__[k + 1 + k * h_dim1], abs(d__1))
1091                                     , d__4 = (d__2 = h__[k + (k + 1) * h_dim1]
1092                                     , abs(d__2));
1093                             h12 = f2cmax(d__3,d__4);
1094 /* Computing MIN */
1095                             d__3 = (d__1 = h__[k + 1 + k * h_dim1], abs(d__1))
1096                                     , d__4 = (d__2 = h__[k + (k + 1) * h_dim1]
1097                                     , abs(d__2));
1098                             h21 = f2cmin(d__3,d__4);
1099 /* Computing MAX */
1100                             d__3 = (d__1 = h__[k + 1 + (k + 1) * h_dim1], abs(
1101                                     d__1)), d__4 = (d__2 = h__[k + k * h_dim1]
1102                                      - h__[k + 1 + (k + 1) * h_dim1], abs(
1103                                     d__2));
1104                             h11 = f2cmax(d__3,d__4);
1105 /* Computing MIN */
1106                             d__3 = (d__1 = h__[k + 1 + (k + 1) * h_dim1], abs(
1107                                     d__1)), d__4 = (d__2 = h__[k + k * h_dim1]
1108                                      - h__[k + 1 + (k + 1) * h_dim1], abs(
1109                                     d__2));
1110                             h22 = f2cmin(d__3,d__4);
1111                             scl = h11 + h12;
1112                             tst2 = h22 * (h11 / scl);
1113
1114 /* Computing MAX */
1115                             d__1 = smlnum, d__2 = ulp * tst2;
1116                             if (tst2 == 0. || h21 * (h12 / scl) <= f2cmax(d__1,
1117                                     d__2)) {
1118                                 h__[k + 1 + k * h_dim1] = 0.;
1119                             }
1120                         }
1121                     }
1122                 }
1123
1124 /*              ==== Accumulate orthogonal transformations. ==== */
1125
1126                 if (accum) {
1127                     kms = k - incol;
1128 /* Computing MAX */
1129                     i__4 = 1, i__5 = *ktop - incol;
1130                     i__6 = kdu;
1131                     for (j = f2cmax(i__4,i__5); j <= i__6; ++j) {
1132                         refsum = v[m22 * v_dim1 + 1] * (u[j + (kms + 1) * 
1133                                 u_dim1] + v[m22 * v_dim1 + 2] * u[j + (kms + 
1134                                 2) * u_dim1]);
1135                         u[j + (kms + 1) * u_dim1] -= refsum;
1136                         u[j + (kms + 2) * u_dim1] -= refsum * v[m22 * v_dim1 
1137                                 + 2];
1138 /* L50: */
1139                     }
1140                 } else if (*wantz) {
1141                     i__6 = *ihiz;
1142                     for (j = *iloz; j <= i__6; ++j) {
1143                         refsum = v[m22 * v_dim1 + 1] * (z__[j + (k + 1) * 
1144                                 z_dim1] + v[m22 * v_dim1 + 2] * z__[j + (k + 
1145                                 2) * z_dim1]);
1146                         z__[j + (k + 1) * z_dim1] -= refsum;
1147                         z__[j + (k + 2) * z_dim1] -= refsum * v[m22 * v_dim1 
1148                                 + 2];
1149 /* L60: */
1150                     }
1151                 }
1152             }
1153
1154 /*           ==== Normal case: Chain of 3-by-3 reflections ==== */
1155
1156             i__6 = mtop;
1157             for (m = mbot; m >= i__6; --m) {
1158                 k = krcol + (m - 1 << 1);
1159                 if (k == *ktop - 1) {
1160                     dlaqr1_(&c__3, &h__[*ktop + *ktop * h_dim1], ldh, &sr[(m 
1161                             << 1) - 1], &si[(m << 1) - 1], &sr[m * 2], &si[m *
1162                              2], &v[m * v_dim1 + 1]);
1163                     alpha = v[m * v_dim1 + 1];
1164                     dlarfg_(&c__3, &alpha, &v[m * v_dim1 + 2], &c__1, &v[m * 
1165                             v_dim1 + 1]);
1166                 } else {
1167
1168 /*                 ==== Perform delayed transformation of row below */
1169 /*                 .    Mth bulge. Exploit fact that first two elements */
1170 /*                 .    of row are actually zero. ==== */
1171
1172                     refsum = v[m * v_dim1 + 1] * v[m * v_dim1 + 3] * h__[k + 
1173                             3 + (k + 2) * h_dim1];
1174                     h__[k + 3 + k * h_dim1] = -refsum;
1175                     h__[k + 3 + (k + 1) * h_dim1] = -refsum * v[m * v_dim1 + 
1176                             2];
1177                     h__[k + 3 + (k + 2) * h_dim1] -= refsum * v[m * v_dim1 + 
1178                             3];
1179
1180 /*                 ==== Calculate reflection to move */
1181 /*                 .    Mth bulge one step. ==== */
1182
1183                     beta = h__[k + 1 + k * h_dim1];
1184                     v[m * v_dim1 + 2] = h__[k + 2 + k * h_dim1];
1185                     v[m * v_dim1 + 3] = h__[k + 3 + k * h_dim1];
1186                     dlarfg_(&c__3, &beta, &v[m * v_dim1 + 2], &c__1, &v[m * 
1187                             v_dim1 + 1]);
1188
1189 /*                 ==== A Bulge may collapse because of vigilant */
1190 /*                 .    deflation or destructive underflow.  In the */
1191 /*                 .    underflow case, try the two-small-subdiagonals */
1192 /*                 .    trick to try to reinflate the bulge.  ==== */
1193
1194                     if (h__[k + 3 + k * h_dim1] != 0. || h__[k + 3 + (k + 1) *
1195                              h_dim1] != 0. || h__[k + 3 + (k + 2) * h_dim1] ==
1196                              0.) {
1197
1198 /*                    ==== Typical case: not collapsed (yet). ==== */
1199
1200                         h__[k + 1 + k * h_dim1] = beta;
1201                         h__[k + 2 + k * h_dim1] = 0.;
1202                         h__[k + 3 + k * h_dim1] = 0.;
1203                     } else {
1204
1205 /*                    ==== Atypical case: collapsed.  Attempt to */
1206 /*                    .    reintroduce ignoring H(K+1,K) and H(K+2,K). */
1207 /*                    .    If the fill resulting from the new */
1208 /*                    .    reflector is too large, then abandon it. */
1209 /*                    .    Otherwise, use the new one. ==== */
1210
1211                         dlaqr1_(&c__3, &h__[k + 1 + (k + 1) * h_dim1], ldh, &
1212                                 sr[(m << 1) - 1], &si[(m << 1) - 1], &sr[m * 
1213                                 2], &si[m * 2], vt);
1214                         alpha = vt[0];
1215                         dlarfg_(&c__3, &alpha, &vt[1], &c__1, vt);
1216                         refsum = vt[0] * (h__[k + 1 + k * h_dim1] + vt[1] * 
1217                                 h__[k + 2 + k * h_dim1]);
1218
1219                         if ((d__1 = h__[k + 2 + k * h_dim1] - refsum * vt[1], 
1220                                 abs(d__1)) + (d__2 = refsum * vt[2], abs(d__2)
1221                                 ) > ulp * ((d__3 = h__[k + k * h_dim1], abs(
1222                                 d__3)) + (d__4 = h__[k + 1 + (k + 1) * h_dim1]
1223                                 , abs(d__4)) + (d__5 = h__[k + 2 + (k + 2) * 
1224                                 h_dim1], abs(d__5)))) {
1225
1226 /*                       ==== Starting a new bulge here would */
1227 /*                       .    create non-negligible fill.  Use */
1228 /*                       .    the old one with trepidation. ==== */
1229
1230                             h__[k + 1 + k * h_dim1] = beta;
1231                             h__[k + 2 + k * h_dim1] = 0.;
1232                             h__[k + 3 + k * h_dim1] = 0.;
1233                         } else {
1234
1235 /*                       ==== Starting a new bulge here would */
1236 /*                       .    create only negligible fill. */
1237 /*                       .    Replace the old reflector with */
1238 /*                       .    the new one. ==== */
1239
1240                             h__[k + 1 + k * h_dim1] -= refsum;
1241                             h__[k + 2 + k * h_dim1] = 0.;
1242                             h__[k + 3 + k * h_dim1] = 0.;
1243                             v[m * v_dim1 + 1] = vt[0];
1244                             v[m * v_dim1 + 2] = vt[1];
1245                             v[m * v_dim1 + 3] = vt[2];
1246                         }
1247                     }
1248                 }
1249
1250 /*              ====  Apply reflection from the right and */
1251 /*              .     the first column of update from the left. */
1252 /*              .     These updates are required for the vigilant */
1253 /*              .     deflation check. We still delay most of the */
1254 /*              .     updates from the left for efficiency. ==== */
1255
1256 /* Computing MIN */
1257                 i__5 = *kbot, i__7 = k + 3;
1258                 i__4 = f2cmin(i__5,i__7);
1259                 for (j = jtop; j <= i__4; ++j) {
1260                     refsum = v[m * v_dim1 + 1] * (h__[j + (k + 1) * h_dim1] + 
1261                             v[m * v_dim1 + 2] * h__[j + (k + 2) * h_dim1] + v[
1262                             m * v_dim1 + 3] * h__[j + (k + 3) * h_dim1]);
1263                     h__[j + (k + 1) * h_dim1] -= refsum;
1264                     h__[j + (k + 2) * h_dim1] -= refsum * v[m * v_dim1 + 2];
1265                     h__[j + (k + 3) * h_dim1] -= refsum * v[m * v_dim1 + 3];
1266 /* L70: */
1267                 }
1268
1269 /*              ==== Perform update from left for subsequent */
1270 /*              .    column. ==== */
1271
1272                 refsum = v[m * v_dim1 + 1] * (h__[k + 1 + (k + 1) * h_dim1] + 
1273                         v[m * v_dim1 + 2] * h__[k + 2 + (k + 1) * h_dim1] + v[
1274                         m * v_dim1 + 3] * h__[k + 3 + (k + 1) * h_dim1]);
1275                 h__[k + 1 + (k + 1) * h_dim1] -= refsum;
1276                 h__[k + 2 + (k + 1) * h_dim1] -= refsum * v[m * v_dim1 + 2];
1277                 h__[k + 3 + (k + 1) * h_dim1] -= refsum * v[m * v_dim1 + 3];
1278
1279 /*              ==== The following convergence test requires that */
1280 /*              .    the tradition small-compared-to-nearby-diagonals */
1281 /*              .    criterion and the Ahues & Tisseur (LAWN 122, 1997) */
1282 /*              .    criteria both be satisfied.  The latter improves */
1283 /*              .    accuracy in some examples. Falling back on an */
1284 /*              .    alternate convergence criterion when TST1 or TST2 */
1285 /*              .    is zero (as done here) is traditional but probably */
1286 /*              .    unnecessary. ==== */
1287
1288                 if (k < *ktop) {
1289                     mycycle_();
1290                 }
1291                 if (h__[k + 1 + k * h_dim1] != 0.) {
1292                     tst1 = (d__1 = h__[k + k * h_dim1], abs(d__1)) + (d__2 = 
1293                             h__[k + 1 + (k + 1) * h_dim1], abs(d__2));
1294                     if (tst1 == 0.) {
1295                         if (k >= *ktop + 1) {
1296                             tst1 += (d__1 = h__[k + (k - 1) * h_dim1], abs(
1297                                     d__1));
1298                         }
1299                         if (k >= *ktop + 2) {
1300                             tst1 += (d__1 = h__[k + (k - 2) * h_dim1], abs(
1301                                     d__1));
1302                         }
1303                         if (k >= *ktop + 3) {
1304                             tst1 += (d__1 = h__[k + (k - 3) * h_dim1], abs(
1305                                     d__1));
1306                         }
1307                         if (k <= *kbot - 2) {
1308                             tst1 += (d__1 = h__[k + 2 + (k + 1) * h_dim1], 
1309                                     abs(d__1));
1310                         }
1311                         if (k <= *kbot - 3) {
1312                             tst1 += (d__1 = h__[k + 3 + (k + 1) * h_dim1], 
1313                                     abs(d__1));
1314                         }
1315                         if (k <= *kbot - 4) {
1316                             tst1 += (d__1 = h__[k + 4 + (k + 1) * h_dim1], 
1317                                     abs(d__1));
1318                         }
1319                     }
1320 /* Computing MAX */
1321                     d__2 = smlnum, d__3 = ulp * tst1;
1322                     if ((d__1 = h__[k + 1 + k * h_dim1], abs(d__1)) <= f2cmax(
1323                             d__2,d__3)) {
1324 /* Computing MAX */
1325                         d__3 = (d__1 = h__[k + 1 + k * h_dim1], abs(d__1)), 
1326                                 d__4 = (d__2 = h__[k + (k + 1) * h_dim1], abs(
1327                                 d__2));
1328                         h12 = f2cmax(d__3,d__4);
1329 /* Computing MIN */
1330                         d__3 = (d__1 = h__[k + 1 + k * h_dim1], abs(d__1)), 
1331                                 d__4 = (d__2 = h__[k + (k + 1) * h_dim1], abs(
1332                                 d__2));
1333                         h21 = f2cmin(d__3,d__4);
1334 /* Computing MAX */
1335                         d__3 = (d__1 = h__[k + 1 + (k + 1) * h_dim1], abs(
1336                                 d__1)), d__4 = (d__2 = h__[k + k * h_dim1] - 
1337                                 h__[k + 1 + (k + 1) * h_dim1], abs(d__2));
1338                         h11 = f2cmax(d__3,d__4);
1339 /* Computing MIN */
1340                         d__3 = (d__1 = h__[k + 1 + (k + 1) * h_dim1], abs(
1341                                 d__1)), d__4 = (d__2 = h__[k + k * h_dim1] - 
1342                                 h__[k + 1 + (k + 1) * h_dim1], abs(d__2));
1343                         h22 = f2cmin(d__3,d__4);
1344                         scl = h11 + h12;
1345                         tst2 = h22 * (h11 / scl);
1346
1347 /* Computing MAX */
1348                         d__1 = smlnum, d__2 = ulp * tst2;
1349                         if (tst2 == 0. || h21 * (h12 / scl) <= f2cmax(d__1,d__2))
1350                                  {
1351                             h__[k + 1 + k * h_dim1] = 0.;
1352                         }
1353                     }
1354                 }
1355 /* L80: */
1356             }
1357
1358 /*           ==== Multiply H by reflections from the left ==== */
1359
1360             if (accum) {
1361                 jbot = f2cmin(ndcol,*kbot);
1362             } else if (*wantt) {
1363                 jbot = *n;
1364             } else {
1365                 jbot = *kbot;
1366             }
1367
1368             i__6 = mtop;
1369             for (m = mbot; m >= i__6; --m) {
1370                 k = krcol + (m - 1 << 1);
1371 /* Computing MAX */
1372                 i__4 = *ktop, i__5 = krcol + (m << 1);
1373                 i__7 = jbot;
1374                 for (j = f2cmax(i__4,i__5); j <= i__7; ++j) {
1375                     refsum = v[m * v_dim1 + 1] * (h__[k + 1 + j * h_dim1] + v[
1376                             m * v_dim1 + 2] * h__[k + 2 + j * h_dim1] + v[m * 
1377                             v_dim1 + 3] * h__[k + 3 + j * h_dim1]);
1378                     h__[k + 1 + j * h_dim1] -= refsum;
1379                     h__[k + 2 + j * h_dim1] -= refsum * v[m * v_dim1 + 2];
1380                     h__[k + 3 + j * h_dim1] -= refsum * v[m * v_dim1 + 3];
1381 /* L90: */
1382                 }
1383 /* L100: */
1384             }
1385
1386 /*           ==== Accumulate orthogonal transformations. ==== */
1387
1388             if (accum) {
1389
1390 /*              ==== Accumulate U. (If needed, update Z later */
1391 /*              .    with an efficient matrix-matrix */
1392 /*              .    multiply.) ==== */
1393
1394                 i__6 = mtop;
1395                 for (m = mbot; m >= i__6; --m) {
1396                     k = krcol + (m - 1 << 1);
1397                     kms = k - incol;
1398 /* Computing MAX */
1399                     i__7 = 1, i__4 = *ktop - incol;
1400                     i2 = f2cmax(i__7,i__4);
1401 /* Computing MAX */
1402                     i__7 = i2, i__4 = kms - (krcol - incol) + 1;
1403                     i2 = f2cmax(i__7,i__4);
1404 /* Computing MIN */
1405                     i__7 = kdu, i__4 = krcol + (mbot - 1 << 1) - incol + 5;
1406                     i4 = f2cmin(i__7,i__4);
1407                     i__7 = i4;
1408                     for (j = i2; j <= i__7; ++j) {
1409                         refsum = v[m * v_dim1 + 1] * (u[j + (kms + 1) * 
1410                                 u_dim1] + v[m * v_dim1 + 2] * u[j + (kms + 2) 
1411                                 * u_dim1] + v[m * v_dim1 + 3] * u[j + (kms + 
1412                                 3) * u_dim1]);
1413                         u[j + (kms + 1) * u_dim1] -= refsum;
1414                         u[j + (kms + 2) * u_dim1] -= refsum * v[m * v_dim1 + 
1415                                 2];
1416                         u[j + (kms + 3) * u_dim1] -= refsum * v[m * v_dim1 + 
1417                                 3];
1418 /* L110: */
1419                     }
1420 /* L120: */
1421                 }
1422             } else if (*wantz) {
1423
1424 /*              ==== U is not accumulated, so update Z */
1425 /*              .    now by multiplying by reflections */
1426 /*              .    from the right. ==== */
1427
1428                 i__6 = mtop;
1429                 for (m = mbot; m >= i__6; --m) {
1430                     k = krcol + (m - 1 << 1);
1431                     i__7 = *ihiz;
1432                     for (j = *iloz; j <= i__7; ++j) {
1433                         refsum = v[m * v_dim1 + 1] * (z__[j + (k + 1) * 
1434                                 z_dim1] + v[m * v_dim1 + 2] * z__[j + (k + 2) 
1435                                 * z_dim1] + v[m * v_dim1 + 3] * z__[j + (k + 
1436                                 3) * z_dim1]);
1437                         z__[j + (k + 1) * z_dim1] -= refsum;
1438                         z__[j + (k + 2) * z_dim1] -= refsum * v[m * v_dim1 + 
1439                                 2];
1440                         z__[j + (k + 3) * z_dim1] -= refsum * v[m * v_dim1 + 
1441                                 3];
1442 /* L130: */
1443                     }
1444 /* L140: */
1445                 }
1446             }
1447
1448 /*           ==== End of near-the-diagonal bulge chase. ==== */
1449
1450 /* L145: */
1451         }
1452
1453 /*        ==== Use U (if accumulated) to update far-from-diagonal */
1454 /*        .    entries in H.  If required, use U to update Z as */
1455 /*        .    well. ==== */
1456
1457         if (accum) {
1458             if (*wantt) {
1459                 jtop = 1;
1460                 jbot = *n;
1461             } else {
1462                 jtop = *ktop;
1463                 jbot = *kbot;
1464             }
1465 /* Computing MAX */
1466             i__3 = 1, i__6 = *ktop - incol;
1467             k1 = f2cmax(i__3,i__6);
1468 /* Computing MAX */
1469             i__3 = 0, i__6 = ndcol - *kbot;
1470             nu = kdu - f2cmax(i__3,i__6) - k1 + 1;
1471
1472 /*           ==== Horizontal Multiply ==== */
1473
1474             i__3 = jbot;
1475             i__6 = *nh;
1476             for (jcol = f2cmin(ndcol,*kbot) + 1; i__6 < 0 ? jcol >= i__3 : jcol 
1477                     <= i__3; jcol += i__6) {
1478 /* Computing MIN */
1479                 i__7 = *nh, i__4 = jbot - jcol + 1;
1480                 jlen = f2cmin(i__7,i__4);
1481                 dgemm_("C", "N", &nu, &jlen, &nu, &c_b8, &u[k1 + k1 * u_dim1],
1482                          ldu, &h__[incol + k1 + jcol * h_dim1], ldh, &c_b7, &
1483                         wh[wh_offset], ldwh);
1484                 dlacpy_("ALL", &nu, &jlen, &wh[wh_offset], ldwh, &h__[incol + 
1485                         k1 + jcol * h_dim1], ldh);
1486 /* L150: */
1487             }
1488
1489 /*           ==== Vertical multiply ==== */
1490
1491             i__6 = f2cmax(*ktop,incol) - 1;
1492             i__3 = *nv;
1493             for (jrow = jtop; i__3 < 0 ? jrow >= i__6 : jrow <= i__6; jrow += 
1494                     i__3) {
1495 /* Computing MIN */
1496                 i__7 = *nv, i__4 = f2cmax(*ktop,incol) - jrow;
1497                 jlen = f2cmin(i__7,i__4);
1498                 dgemm_("N", "N", &jlen, &nu, &nu, &c_b8, &h__[jrow + (incol + 
1499                         k1) * h_dim1], ldh, &u[k1 + k1 * u_dim1], ldu, &c_b7, 
1500                         &wv[wv_offset], ldwv);
1501                 dlacpy_("ALL", &jlen, &nu, &wv[wv_offset], ldwv, &h__[jrow + (
1502                         incol + k1) * h_dim1], ldh);
1503 /* L160: */
1504             }
1505
1506 /*           ==== Z multiply (also vertical) ==== */
1507
1508             if (*wantz) {
1509                 i__3 = *ihiz;
1510                 i__6 = *nv;
1511                 for (jrow = *iloz; i__6 < 0 ? jrow >= i__3 : jrow <= i__3; 
1512                         jrow += i__6) {
1513 /* Computing MIN */
1514                     i__7 = *nv, i__4 = *ihiz - jrow + 1;
1515                     jlen = f2cmin(i__7,i__4);
1516                     dgemm_("N", "N", &jlen, &nu, &nu, &c_b8, &z__[jrow + (
1517                             incol + k1) * z_dim1], ldz, &u[k1 + k1 * u_dim1], 
1518                             ldu, &c_b7, &wv[wv_offset], ldwv);
1519                     dlacpy_("ALL", &jlen, &nu, &wv[wv_offset], ldwv, &z__[
1520                             jrow + (incol + k1) * z_dim1], ldz);
1521 /* L170: */
1522                 }
1523             }
1524         }
1525 /* L180: */
1526     }
1527
1528 /*     ==== End of DLAQR5 ==== */
1529
1530     return 0;
1531 } /* dlaqr5_ */
1532