C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / claqr2.c
1 #include <math.h>
2 #include <stdlib.h>
3 #include <string.h>
4 #include <stdio.h>
5 #include <complex.h>
6 #ifdef complex
7 #undef complex
8 #endif
9 #ifdef I
10 #undef I
11 #endif
12
13 #if defined(_WIN64)
14 typedef long long BLASLONG;
15 typedef unsigned long long BLASULONG;
16 #else
17 typedef long BLASLONG;
18 typedef unsigned long BLASULONG;
19 #endif
20
21 #ifdef LAPACK_ILP64
22 typedef BLASLONG blasint;
23 #if defined(_WIN64)
24 #define blasabs(x) llabs(x)
25 #else
26 #define blasabs(x) labs(x)
27 #endif
28 #else
29 typedef int blasint;
30 #define blasabs(x) abs(x)
31 #endif
32
33 typedef blasint integer;
34
35 typedef unsigned int uinteger;
36 typedef char *address;
37 typedef short int shortint;
38 typedef float real;
39 typedef double doublereal;
40 typedef struct { real r, i; } complex;
41 typedef struct { doublereal r, i; } doublecomplex;
42 #ifdef _MSC_VER
43 static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
44 static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
45 static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
46 static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
47 #else
48 static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
49 static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
50 static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
51 static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
52 #endif
53 #define pCf(z) (*_pCf(z))
54 #define pCd(z) (*_pCd(z))
55 typedef int logical;
56 typedef short int shortlogical;
57 typedef char logical1;
58 typedef char integer1;
59
60 #define TRUE_ (1)
61 #define FALSE_ (0)
62
63 /* Extern is for use with -E */
64 #ifndef Extern
65 #define Extern extern
66 #endif
67
68 /* I/O stuff */
69
70 typedef int flag;
71 typedef int ftnlen;
72 typedef int ftnint;
73
74 /*external read, write*/
75 typedef struct
76 {       flag cierr;
77         ftnint ciunit;
78         flag ciend;
79         char *cifmt;
80         ftnint cirec;
81 } cilist;
82
83 /*internal read, write*/
84 typedef struct
85 {       flag icierr;
86         char *iciunit;
87         flag iciend;
88         char *icifmt;
89         ftnint icirlen;
90         ftnint icirnum;
91 } icilist;
92
93 /*open*/
94 typedef struct
95 {       flag oerr;
96         ftnint ounit;
97         char *ofnm;
98         ftnlen ofnmlen;
99         char *osta;
100         char *oacc;
101         char *ofm;
102         ftnint orl;
103         char *oblnk;
104 } olist;
105
106 /*close*/
107 typedef struct
108 {       flag cerr;
109         ftnint cunit;
110         char *csta;
111 } cllist;
112
113 /*rewind, backspace, endfile*/
114 typedef struct
115 {       flag aerr;
116         ftnint aunit;
117 } alist;
118
119 /* inquire */
120 typedef struct
121 {       flag inerr;
122         ftnint inunit;
123         char *infile;
124         ftnlen infilen;
125         ftnint  *inex;  /*parameters in standard's order*/
126         ftnint  *inopen;
127         ftnint  *innum;
128         ftnint  *innamed;
129         char    *inname;
130         ftnlen  innamlen;
131         char    *inacc;
132         ftnlen  inacclen;
133         char    *inseq;
134         ftnlen  inseqlen;
135         char    *indir;
136         ftnlen  indirlen;
137         char    *infmt;
138         ftnlen  infmtlen;
139         char    *inform;
140         ftnint  informlen;
141         char    *inunf;
142         ftnlen  inunflen;
143         ftnint  *inrecl;
144         ftnint  *innrec;
145         char    *inblank;
146         ftnlen  inblanklen;
147 } inlist;
148
149 #define VOID void
150
151 union Multitype {       /* for multiple entry points */
152         integer1 g;
153         shortint h;
154         integer i;
155         /* longint j; */
156         real r;
157         doublereal d;
158         complex c;
159         doublecomplex z;
160         };
161
162 typedef union Multitype Multitype;
163
164 struct Vardesc {        /* for Namelist */
165         char *name;
166         char *addr;
167         ftnlen *dims;
168         int  type;
169         };
170 typedef struct Vardesc Vardesc;
171
172 struct Namelist {
173         char *name;
174         Vardesc **vars;
175         int nvars;
176         };
177 typedef struct Namelist Namelist;
178
179 #define abs(x) ((x) >= 0 ? (x) : -(x))
180 #define dabs(x) (fabs(x))
181 #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
182 #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
183 #define dmin(a,b) (f2cmin(a,b))
184 #define dmax(a,b) (f2cmax(a,b))
185 #define bit_test(a,b)   ((a) >> (b) & 1)
186 #define bit_clear(a,b)  ((a) & ~((uinteger)1 << (b)))
187 #define bit_set(a,b)    ((a) |  ((uinteger)1 << (b)))
188
189 #define abort_() { sig_die("Fortran abort routine called", 1); }
190 #define c_abs(z) (cabsf(Cf(z)))
191 #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
192 #ifdef _MSC_VER
193 #define c_div(c, a, b) {Cf(c)._Val[0] = (Cf(a)._Val[0]/Cf(b)._Val[0]); Cf(c)._Val[1]=(Cf(a)._Val[1]/Cf(b)._Val[1]);}
194 #define z_div(c, a, b) {Cd(c)._Val[0] = (Cd(a)._Val[0]/Cd(b)._Val[0]); Cd(c)._Val[1]=(Cd(a)._Val[1]/df(b)._Val[1]);}
195 #else
196 #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
197 #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
198 #endif
199 #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
200 #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
201 #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
202 //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
203 #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
204 #define d_abs(x) (fabs(*(x)))
205 #define d_acos(x) (acos(*(x)))
206 #define d_asin(x) (asin(*(x)))
207 #define d_atan(x) (atan(*(x)))
208 #define d_atn2(x, y) (atan2(*(x),*(y)))
209 #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
210 #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
211 #define d_cos(x) (cos(*(x)))
212 #define d_cosh(x) (cosh(*(x)))
213 #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
214 #define d_exp(x) (exp(*(x)))
215 #define d_imag(z) (cimag(Cd(z)))
216 #define r_imag(z) (cimagf(Cf(z)))
217 #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
218 #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
219 #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
220 #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
221 #define d_log(x) (log(*(x)))
222 #define d_mod(x, y) (fmod(*(x), *(y)))
223 #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
224 #define d_nint(x) u_nint(*(x))
225 #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
226 #define d_sign(a,b) u_sign(*(a),*(b))
227 #define r_sign(a,b) u_sign(*(a),*(b))
228 #define d_sin(x) (sin(*(x)))
229 #define d_sinh(x) (sinh(*(x)))
230 #define d_sqrt(x) (sqrt(*(x)))
231 #define d_tan(x) (tan(*(x)))
232 #define d_tanh(x) (tanh(*(x)))
233 #define i_abs(x) abs(*(x))
234 #define i_dnnt(x) ((integer)u_nint(*(x)))
235 #define i_len(s, n) (n)
236 #define i_nint(x) ((integer)u_nint(*(x)))
237 #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
238 #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
239 #define pow_si(B,E) spow_ui(*(B),*(E))
240 #define pow_ri(B,E) spow_ui(*(B),*(E))
241 #define pow_di(B,E) dpow_ui(*(B),*(E))
242 #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
243 #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
244 #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
245 #define s_cat(lpp, rpp, rnp, np, llp) {         ftnlen i, nc, ll; char *f__rp, *lp;     ll = (llp); lp = (lpp);         for(i=0; i < (int)*(np); ++i) {                 nc = ll;                if((rnp)[i] < nc) nc = (rnp)[i];                ll -= nc;               f__rp = (rpp)[i];               while(--nc >= 0) *lp++ = *(f__rp)++;         }  while(--ll >= 0) *lp++ = ' '; }
246 #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
247 #define s_copy(A,B,C,D) { int __i,__m; for (__i=0, __m=f2cmin((C),(D)); __i<__m && (B)[__i] != 0; ++__i) (A)[__i] = (B)[__i]; }
248 #define sig_die(s, kill) { exit(1); }
249 #define s_stop(s, n) {exit(0);}
250 static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
251 #define z_abs(z) (cabs(Cd(z)))
252 #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
253 #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
254 #define myexit_() break;
255 #define mycycle() continue;
256 #define myceiling(w) {ceil(w)}
257 #define myhuge(w) {HUGE_VAL}
258 //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
259 #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
260
261 /* procedure parameter types for -A and -C++ */
262
263 #define F2C_proc_par_types 1
264 #ifdef __cplusplus
265 typedef logical (*L_fp)(...);
266 #else
267 typedef logical (*L_fp)();
268 #endif
269
270 static float spow_ui(float x, integer n) {
271         float pow=1.0; unsigned long int u;
272         if(n != 0) {
273                 if(n < 0) n = -n, x = 1/x;
274                 for(u = n; ; ) {
275                         if(u & 01) pow *= x;
276                         if(u >>= 1) x *= x;
277                         else break;
278                 }
279         }
280         return pow;
281 }
282 static double dpow_ui(double x, integer n) {
283         double pow=1.0; unsigned long int u;
284         if(n != 0) {
285                 if(n < 0) n = -n, x = 1/x;
286                 for(u = n; ; ) {
287                         if(u & 01) pow *= x;
288                         if(u >>= 1) x *= x;
289                         else break;
290                 }
291         }
292         return pow;
293 }
294 #ifdef _MSC_VER
295 static _Fcomplex cpow_ui(complex x, integer n) {
296         complex pow={1.0,0.0}; unsigned long int u;
297                 if(n != 0) {
298                 if(n < 0) n = -n, x.r = 1/x.r, x.i=1/x.i;
299                 for(u = n; ; ) {
300                         if(u & 01) pow.r *= x.r, pow.i *= x.i;
301                         if(u >>= 1) x.r *= x.r, x.i *= x.i;
302                         else break;
303                 }
304         }
305         _Fcomplex p={pow.r, pow.i};
306         return p;
307 }
308 #else
309 static _Complex float cpow_ui(_Complex float x, integer n) {
310         _Complex float pow=1.0; unsigned long int u;
311         if(n != 0) {
312                 if(n < 0) n = -n, x = 1/x;
313                 for(u = n; ; ) {
314                         if(u & 01) pow *= x;
315                         if(u >>= 1) x *= x;
316                         else break;
317                 }
318         }
319         return pow;
320 }
321 #endif
322 #ifdef _MSC_VER
323 static _Dcomplex zpow_ui(_Dcomplex x, integer n) {
324         _Dcomplex pow={1.0,0.0}; unsigned long int u;
325         if(n != 0) {
326                 if(n < 0) n = -n, x._Val[0] = 1/x._Val[0], x._Val[1] =1/x._Val[1];
327                 for(u = n; ; ) {
328                         if(u & 01) pow._Val[0] *= x._Val[0], pow._Val[1] *= x._Val[1];
329                         if(u >>= 1) x._Val[0] *= x._Val[0], x._Val[1] *= x._Val[1];
330                         else break;
331                 }
332         }
333         _Dcomplex p = {pow._Val[0], pow._Val[1]};
334         return p;
335 }
336 #else
337 static _Complex double zpow_ui(_Complex double x, integer n) {
338         _Complex double pow=1.0; unsigned long int u;
339         if(n != 0) {
340                 if(n < 0) n = -n, x = 1/x;
341                 for(u = n; ; ) {
342                         if(u & 01) pow *= x;
343                         if(u >>= 1) x *= x;
344                         else break;
345                 }
346         }
347         return pow;
348 }
349 #endif
350 static integer pow_ii(integer x, integer n) {
351         integer pow; unsigned long int u;
352         if (n <= 0) {
353                 if (n == 0 || x == 1) pow = 1;
354                 else if (x != -1) pow = x == 0 ? 1/x : 0;
355                 else n = -n;
356         }
357         if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
358                 u = n;
359                 for(pow = 1; ; ) {
360                         if(u & 01) pow *= x;
361                         if(u >>= 1) x *= x;
362                         else break;
363                 }
364         }
365         return pow;
366 }
367 static integer dmaxloc_(double *w, integer s, integer e, integer *n)
368 {
369         double m; integer i, mi;
370         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
371                 if (w[i-1]>m) mi=i ,m=w[i-1];
372         return mi-s+1;
373 }
374 static integer smaxloc_(float *w, integer s, integer e, integer *n)
375 {
376         float m; integer i, mi;
377         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
378                 if (w[i-1]>m) mi=i ,m=w[i-1];
379         return mi-s+1;
380 }
381 static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
382         integer n = *n_, incx = *incx_, incy = *incy_, i;
383 #ifdef _MSC_VER
384         _Fcomplex zdotc = {0.0, 0.0};
385         if (incx == 1 && incy == 1) {
386                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
387                         zdotc._Val[0] += conjf(Cf(&x[i]))._Val[0] * Cf(&y[i])._Val[0];
388                         zdotc._Val[1] += conjf(Cf(&x[i]))._Val[1] * Cf(&y[i])._Val[1];
389                 }
390         } else {
391                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
392                         zdotc._Val[0] += conjf(Cf(&x[i*incx]))._Val[0] * Cf(&y[i*incy])._Val[0];
393                         zdotc._Val[1] += conjf(Cf(&x[i*incx]))._Val[1] * Cf(&y[i*incy])._Val[1];
394                 }
395         }
396         pCf(z) = zdotc;
397 }
398 #else
399         _Complex float zdotc = 0.0;
400         if (incx == 1 && incy == 1) {
401                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
402                         zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
403                 }
404         } else {
405                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
406                         zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
407                 }
408         }
409         pCf(z) = zdotc;
410 }
411 #endif
412 static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
413         integer n = *n_, incx = *incx_, incy = *incy_, i;
414 #ifdef _MSC_VER
415         _Dcomplex zdotc = {0.0, 0.0};
416         if (incx == 1 && incy == 1) {
417                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
418                         zdotc._Val[0] += conj(Cd(&x[i]))._Val[0] * Cd(&y[i])._Val[0];
419                         zdotc._Val[1] += conj(Cd(&x[i]))._Val[1] * Cd(&y[i])._Val[1];
420                 }
421         } else {
422                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
423                         zdotc._Val[0] += conj(Cd(&x[i*incx]))._Val[0] * Cd(&y[i*incy])._Val[0];
424                         zdotc._Val[1] += conj(Cd(&x[i*incx]))._Val[1] * Cd(&y[i*incy])._Val[1];
425                 }
426         }
427         pCd(z) = zdotc;
428 }
429 #else
430         _Complex double zdotc = 0.0;
431         if (incx == 1 && incy == 1) {
432                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
433                         zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
434                 }
435         } else {
436                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
437                         zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
438                 }
439         }
440         pCd(z) = zdotc;
441 }
442 #endif  
443 static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
444         integer n = *n_, incx = *incx_, incy = *incy_, i;
445 #ifdef _MSC_VER
446         _Fcomplex zdotc = {0.0, 0.0};
447         if (incx == 1 && incy == 1) {
448                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
449                         zdotc._Val[0] += Cf(&x[i])._Val[0] * Cf(&y[i])._Val[0];
450                         zdotc._Val[1] += Cf(&x[i])._Val[1] * Cf(&y[i])._Val[1];
451                 }
452         } else {
453                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
454                         zdotc._Val[0] += Cf(&x[i*incx])._Val[0] * Cf(&y[i*incy])._Val[0];
455                         zdotc._Val[1] += Cf(&x[i*incx])._Val[1] * Cf(&y[i*incy])._Val[1];
456                 }
457         }
458         pCf(z) = zdotc;
459 }
460 #else
461         _Complex float zdotc = 0.0;
462         if (incx == 1 && incy == 1) {
463                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
464                         zdotc += Cf(&x[i]) * Cf(&y[i]);
465                 }
466         } else {
467                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
468                         zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
469                 }
470         }
471         pCf(z) = zdotc;
472 }
473 #endif
474 static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
475         integer n = *n_, incx = *incx_, incy = *incy_, i;
476 #ifdef _MSC_VER
477         _Dcomplex zdotc = {0.0, 0.0};
478         if (incx == 1 && incy == 1) {
479                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
480                         zdotc._Val[0] += Cd(&x[i])._Val[0] * Cd(&y[i])._Val[0];
481                         zdotc._Val[1] += Cd(&x[i])._Val[1] * Cd(&y[i])._Val[1];
482                 }
483         } else {
484                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
485                         zdotc._Val[0] += Cd(&x[i*incx])._Val[0] * Cd(&y[i*incy])._Val[0];
486                         zdotc._Val[1] += Cd(&x[i*incx])._Val[1] * Cd(&y[i*incy])._Val[1];
487                 }
488         }
489         pCd(z) = zdotc;
490 }
491 #else
492         _Complex double zdotc = 0.0;
493         if (incx == 1 && incy == 1) {
494                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
495                         zdotc += Cd(&x[i]) * Cd(&y[i]);
496                 }
497         } else {
498                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
499                         zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
500                 }
501         }
502         pCd(z) = zdotc;
503 }
504 #endif
505 /*  -- translated by f2c (version 20000121).
506    You must link the resulting object file with the libraries:
507         -lf2c -lm   (in that order)
508 */
509
510
511
512
513 /* Table of constant values */
514
515 static complex c_b1 = {0.f,0.f};
516 static complex c_b2 = {1.f,0.f};
517 static integer c__1 = 1;
518 static integer c_n1 = -1;
519 static logical c_true = TRUE_;
520
521 /* > \brief \b CLAQR2 performs the unitary similarity transformation of a Hessenberg matrix to detect and defl
522 ate fully converged eigenvalues from a trailing principal submatrix (aggressive early deflation). */
523
524 /*  =========== DOCUMENTATION =========== */
525
526 /* Online html documentation available at */
527 /*            http://www.netlib.org/lapack/explore-html/ */
528
529 /* > \htmlonly */
530 /* > Download CLAQR2 + dependencies */
531 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/claqr2.
532 f"> */
533 /* > [TGZ]</a> */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/claqr2.
535 f"> */
536 /* > [ZIP]</a> */
537 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/claqr2.
538 f"> */
539 /* > [TXT]</a> */
540 /* > \endhtmlonly */
541
542 /*  Definition: */
543 /*  =========== */
544
545 /*       SUBROUTINE CLAQR2( WANTT, WANTZ, N, KTOP, KBOT, NW, H, LDH, ILOZ, */
546 /*                          IHIZ, Z, LDZ, NS, ND, SH, V, LDV, NH, T, LDT, */
547 /*                          NV, WV, LDWV, WORK, LWORK ) */
548
549 /*       INTEGER            IHIZ, ILOZ, KBOT, KTOP, LDH, LDT, LDV, LDWV, */
550 /*      $                   LDZ, LWORK, N, ND, NH, NS, NV, NW */
551 /*       LOGICAL            WANTT, WANTZ */
552 /*       COMPLEX            H( LDH, * ), SH( * ), T( LDT, * ), V( LDV, * ), */
553 /*      $                   WORK( * ), WV( LDWV, * ), Z( LDZ, * ) */
554
555
556 /* > \par Purpose: */
557 /*  ============= */
558 /* > */
559 /* > \verbatim */
560 /* > */
561 /* >    CLAQR2 is identical to CLAQR3 except that it avoids */
562 /* >    recursion by calling CLAHQR instead of CLAQR4. */
563 /* > */
564 /* >    Aggressive early deflation: */
565 /* > */
566 /* >    This subroutine accepts as input an upper Hessenberg matrix */
567 /* >    H and performs an unitary similarity transformation */
568 /* >    designed to detect and deflate fully converged eigenvalues from */
569 /* >    a trailing principal submatrix.  On output H has been over- */
570 /* >    written by a new Hessenberg matrix that is a perturbation of */
571 /* >    an unitary similarity transformation of H.  It is to be */
572 /* >    hoped that the final version of H has many zero subdiagonal */
573 /* >    entries. */
574 /* > \endverbatim */
575
576 /*  Arguments: */
577 /*  ========== */
578
579 /* > \param[in] WANTT */
580 /* > \verbatim */
581 /* >          WANTT is LOGICAL */
582 /* >          If .TRUE., then the Hessenberg matrix H is fully updated */
583 /* >          so that the triangular Schur factor may be */
584 /* >          computed (in cooperation with the calling subroutine). */
585 /* >          If .FALSE., then only enough of H is updated to preserve */
586 /* >          the eigenvalues. */
587 /* > \endverbatim */
588 /* > */
589 /* > \param[in] WANTZ */
590 /* > \verbatim */
591 /* >          WANTZ is LOGICAL */
592 /* >          If .TRUE., then the unitary matrix Z is updated so */
593 /* >          so that the unitary Schur factor may be computed */
594 /* >          (in cooperation with the calling subroutine). */
595 /* >          If .FALSE., then Z is not referenced. */
596 /* > \endverbatim */
597 /* > */
598 /* > \param[in] N */
599 /* > \verbatim */
600 /* >          N is INTEGER */
601 /* >          The order of the matrix H and (if WANTZ is .TRUE.) the */
602 /* >          order of the unitary matrix Z. */
603 /* > \endverbatim */
604 /* > */
605 /* > \param[in] KTOP */
606 /* > \verbatim */
607 /* >          KTOP is INTEGER */
608 /* >          It is assumed that either KTOP = 1 or H(KTOP,KTOP-1)=0. */
609 /* >          KBOT and KTOP together determine an isolated block */
610 /* >          along the diagonal of the Hessenberg matrix. */
611 /* > \endverbatim */
612 /* > */
613 /* > \param[in] KBOT */
614 /* > \verbatim */
615 /* >          KBOT is INTEGER */
616 /* >          It is assumed without a check that either */
617 /* >          KBOT = N or H(KBOT+1,KBOT)=0.  KBOT and KTOP together */
618 /* >          determine an isolated block along the diagonal of the */
619 /* >          Hessenberg matrix. */
620 /* > \endverbatim */
621 /* > */
622 /* > \param[in] NW */
623 /* > \verbatim */
624 /* >          NW is INTEGER */
625 /* >          Deflation window size.  1 <= NW <= (KBOT-KTOP+1). */
626 /* > \endverbatim */
627 /* > */
628 /* > \param[in,out] H */
629 /* > \verbatim */
630 /* >          H is COMPLEX array, dimension (LDH,N) */
631 /* >          On input the initial N-by-N section of H stores the */
632 /* >          Hessenberg matrix undergoing aggressive early deflation. */
633 /* >          On output H has been transformed by a unitary */
634 /* >          similarity transformation, perturbed, and the returned */
635 /* >          to Hessenberg form that (it is to be hoped) has some */
636 /* >          zero subdiagonal entries. */
637 /* > \endverbatim */
638 /* > */
639 /* > \param[in] LDH */
640 /* > \verbatim */
641 /* >          LDH is INTEGER */
642 /* >          Leading dimension of H just as declared in the calling */
643 /* >          subroutine.  N <= LDH */
644 /* > \endverbatim */
645 /* > */
646 /* > \param[in] ILOZ */
647 /* > \verbatim */
648 /* >          ILOZ is INTEGER */
649 /* > \endverbatim */
650 /* > */
651 /* > \param[in] IHIZ */
652 /* > \verbatim */
653 /* >          IHIZ is INTEGER */
654 /* >          Specify the rows of Z to which transformations must be */
655 /* >          applied if WANTZ is .TRUE.. 1 <= ILOZ <= IHIZ <= N. */
656 /* > \endverbatim */
657 /* > */
658 /* > \param[in,out] Z */
659 /* > \verbatim */
660 /* >          Z is COMPLEX array, dimension (LDZ,N) */
661 /* >          IF WANTZ is .TRUE., then on output, the unitary */
662 /* >          similarity transformation mentioned above has been */
663 /* >          accumulated into Z(ILOZ:IHIZ,ILOZ:IHIZ) from the right. */
664 /* >          If WANTZ is .FALSE., then Z is unreferenced. */
665 /* > \endverbatim */
666 /* > */
667 /* > \param[in] LDZ */
668 /* > \verbatim */
669 /* >          LDZ is INTEGER */
670 /* >          The leading dimension of Z just as declared in the */
671 /* >          calling subroutine.  1 <= LDZ. */
672 /* > \endverbatim */
673 /* > */
674 /* > \param[out] NS */
675 /* > \verbatim */
676 /* >          NS is INTEGER */
677 /* >          The number of unconverged (ie approximate) eigenvalues */
678 /* >          returned in SR and SI that may be used as shifts by the */
679 /* >          calling subroutine. */
680 /* > \endverbatim */
681 /* > */
682 /* > \param[out] ND */
683 /* > \verbatim */
684 /* >          ND is INTEGER */
685 /* >          The number of converged eigenvalues uncovered by this */
686 /* >          subroutine. */
687 /* > \endverbatim */
688 /* > */
689 /* > \param[out] SH */
690 /* > \verbatim */
691 /* >          SH is COMPLEX array, dimension (KBOT) */
692 /* >          On output, approximate eigenvalues that may */
693 /* >          be used for shifts are stored in SH(KBOT-ND-NS+1) */
694 /* >          through SR(KBOT-ND).  Converged eigenvalues are */
695 /* >          stored in SH(KBOT-ND+1) through SH(KBOT). */
696 /* > \endverbatim */
697 /* > */
698 /* > \param[out] V */
699 /* > \verbatim */
700 /* >          V is COMPLEX array, dimension (LDV,NW) */
701 /* >          An NW-by-NW work array. */
702 /* > \endverbatim */
703 /* > */
704 /* > \param[in] LDV */
705 /* > \verbatim */
706 /* >          LDV is INTEGER */
707 /* >          The leading dimension of V just as declared in the */
708 /* >          calling subroutine.  NW <= LDV */
709 /* > \endverbatim */
710 /* > */
711 /* > \param[in] NH */
712 /* > \verbatim */
713 /* >          NH is INTEGER */
714 /* >          The number of columns of T.  NH >= NW. */
715 /* > \endverbatim */
716 /* > */
717 /* > \param[out] T */
718 /* > \verbatim */
719 /* >          T is COMPLEX array, dimension (LDT,NW) */
720 /* > \endverbatim */
721 /* > */
722 /* > \param[in] LDT */
723 /* > \verbatim */
724 /* >          LDT is INTEGER */
725 /* >          The leading dimension of T just as declared in the */
726 /* >          calling subroutine.  NW <= LDT */
727 /* > \endverbatim */
728 /* > */
729 /* > \param[in] NV */
730 /* > \verbatim */
731 /* >          NV is INTEGER */
732 /* >          The number of rows of work array WV available for */
733 /* >          workspace.  NV >= NW. */
734 /* > \endverbatim */
735 /* > */
736 /* > \param[out] WV */
737 /* > \verbatim */
738 /* >          WV is COMPLEX array, dimension (LDWV,NW) */
739 /* > \endverbatim */
740 /* > */
741 /* > \param[in] LDWV */
742 /* > \verbatim */
743 /* >          LDWV is INTEGER */
744 /* >          The leading dimension of W just as declared in the */
745 /* >          calling subroutine.  NW <= LDV */
746 /* > \endverbatim */
747 /* > */
748 /* > \param[out] WORK */
749 /* > \verbatim */
750 /* >          WORK is COMPLEX array, dimension (LWORK) */
751 /* >          On exit, WORK(1) is set to an estimate of the optimal value */
752 /* >          of LWORK for the given values of N, NW, KTOP and KBOT. */
753 /* > \endverbatim */
754 /* > */
755 /* > \param[in] LWORK */
756 /* > \verbatim */
757 /* >          LWORK is INTEGER */
758 /* >          The dimension of the work array WORK.  LWORK = 2*NW */
759 /* >          suffices, but greater efficiency may result from larger */
760 /* >          values of LWORK. */
761 /* > */
762 /* >          If LWORK = -1, then a workspace query is assumed; CLAQR2 */
763 /* >          only estimates the optimal workspace size for the given */
764 /* >          values of N, NW, KTOP and KBOT.  The estimate is returned */
765 /* >          in WORK(1).  No error message related to LWORK is issued */
766 /* >          by XERBLA.  Neither H nor Z are accessed. */
767 /* > \endverbatim */
768
769 /*  Authors: */
770 /*  ======== */
771
772 /* > \author Univ. of Tennessee */
773 /* > \author Univ. of California Berkeley */
774 /* > \author Univ. of Colorado Denver */
775 /* > \author NAG Ltd. */
776
777 /* > \date June 2017 */
778
779 /* > \ingroup complexOTHERauxiliary */
780
781 /* > \par Contributors: */
782 /*  ================== */
783 /* > */
784 /* >       Karen Braman and Ralph Byers, Department of Mathematics, */
785 /* >       University of Kansas, USA */
786 /* > */
787 /*  ===================================================================== */
788 /* Subroutine */ int claqr2_(logical *wantt, logical *wantz, integer *n, 
789         integer *ktop, integer *kbot, integer *nw, complex *h__, integer *ldh,
790          integer *iloz, integer *ihiz, complex *z__, integer *ldz, integer *
791         ns, integer *nd, complex *sh, complex *v, integer *ldv, integer *nh, 
792         complex *t, integer *ldt, integer *nv, complex *wv, integer *ldwv, 
793         complex *work, integer *lwork)
794 {
795     /* System generated locals */
796     integer h_dim1, h_offset, t_dim1, t_offset, v_dim1, v_offset, wv_dim1, 
797             wv_offset, z_dim1, z_offset, i__1, i__2, i__3, i__4;
798     real r__1, r__2, r__3, r__4, r__5, r__6;
799     complex q__1, q__2;
800
801     /* Local variables */
802     complex beta;
803     integer kcol, info, ifst, ilst, ltop, krow, i__, j;
804     complex s;
805     extern /* Subroutine */ int clarf_(char *, integer *, integer *, complex *
806             , integer *, complex *, complex *, integer *, complex *), 
807             cgemm_(char *, char *, integer *, integer *, integer *, complex *,
808              complex *, integer *, complex *, integer *, complex *, complex *,
809              integer *), ccopy_(integer *, complex *, integer 
810             *, complex *, integer *);
811     integer infqr, kwtop;
812     extern /* Subroutine */ int slabad_(real *, real *), cgehrd_(integer *, 
813             integer *, integer *, complex *, integer *, complex *, complex *, 
814             integer *, integer *), clarfg_(integer *, complex *, complex *, 
815             integer *, complex *);
816     integer jw;
817     extern real slamch_(char *);
818     extern /* Subroutine */ int clahqr_(logical *, logical *, integer *, 
819             integer *, integer *, complex *, integer *, complex *, integer *, 
820             integer *, complex *, integer *, integer *), clacpy_(char *, 
821             integer *, integer *, complex *, integer *, complex *, integer *), claset_(char *, integer *, integer *, complex *, complex 
822             *, complex *, integer *);
823     real safmin, safmax;
824     extern /* Subroutine */ int ctrexc_(char *, integer *, complex *, integer 
825             *, complex *, integer *, integer *, integer *, integer *),
826              cunmhr_(char *, char *, integer *, integer *, integer *, integer 
827             *, complex *, integer *, complex *, complex *, integer *, complex 
828             *, integer *, integer *);
829     real smlnum;
830     integer lwkopt;
831     real foo;
832     integer kln;
833     complex tau;
834     integer knt;
835     real ulp;
836     integer lwk1, lwk2;
837
838
839 /*  -- LAPACK auxiliary routine (version 3.7.1) -- */
840 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
841 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
842 /*     June 2017 */
843
844
845 /*  ================================================================ */
846
847
848 /*     ==== Estimate optimal workspace. ==== */
849
850     /* Parameter adjustments */
851     h_dim1 = *ldh;
852     h_offset = 1 + h_dim1 * 1;
853     h__ -= h_offset;
854     z_dim1 = *ldz;
855     z_offset = 1 + z_dim1 * 1;
856     z__ -= z_offset;
857     --sh;
858     v_dim1 = *ldv;
859     v_offset = 1 + v_dim1 * 1;
860     v -= v_offset;
861     t_dim1 = *ldt;
862     t_offset = 1 + t_dim1 * 1;
863     t -= t_offset;
864     wv_dim1 = *ldwv;
865     wv_offset = 1 + wv_dim1 * 1;
866     wv -= wv_offset;
867     --work;
868
869     /* Function Body */
870 /* Computing MIN */
871     i__1 = *nw, i__2 = *kbot - *ktop + 1;
872     jw = f2cmin(i__1,i__2);
873     if (jw <= 2) {
874         lwkopt = 1;
875     } else {
876
877 /*        ==== Workspace query call to CGEHRD ==== */
878
879         i__1 = jw - 1;
880         cgehrd_(&jw, &c__1, &i__1, &t[t_offset], ldt, &work[1], &work[1], &
881                 c_n1, &info);
882         lwk1 = (integer) work[1].r;
883
884 /*        ==== Workspace query call to CUNMHR ==== */
885
886         i__1 = jw - 1;
887         cunmhr_("R", "N", &jw, &jw, &c__1, &i__1, &t[t_offset], ldt, &work[1],
888                  &v[v_offset], ldv, &work[1], &c_n1, &info);
889         lwk2 = (integer) work[1].r;
890
891 /*        ==== Optimal workspace ==== */
892
893         lwkopt = jw + f2cmax(lwk1,lwk2);
894     }
895
896 /*     ==== Quick return in case of workspace query. ==== */
897
898     if (*lwork == -1) {
899         r__1 = (real) lwkopt;
900         q__1.r = r__1, q__1.i = 0.f;
901         work[1].r = q__1.r, work[1].i = q__1.i;
902         return 0;
903     }
904
905 /*     ==== Nothing to do ... */
906 /*     ... for an empty active block ... ==== */
907     *ns = 0;
908     *nd = 0;
909     work[1].r = 1.f, work[1].i = 0.f;
910     if (*ktop > *kbot) {
911         return 0;
912     }
913 /*     ... nor for an empty deflation window. ==== */
914     if (*nw < 1) {
915         return 0;
916     }
917
918 /*     ==== Machine constants ==== */
919
920     safmin = slamch_("SAFE MINIMUM");
921     safmax = 1.f / safmin;
922     slabad_(&safmin, &safmax);
923     ulp = slamch_("PRECISION");
924     smlnum = safmin * ((real) (*n) / ulp);
925
926 /*     ==== Setup deflation window ==== */
927
928 /* Computing MIN */
929     i__1 = *nw, i__2 = *kbot - *ktop + 1;
930     jw = f2cmin(i__1,i__2);
931     kwtop = *kbot - jw + 1;
932     if (kwtop == *ktop) {
933         s.r = 0.f, s.i = 0.f;
934     } else {
935         i__1 = kwtop + (kwtop - 1) * h_dim1;
936         s.r = h__[i__1].r, s.i = h__[i__1].i;
937     }
938
939     if (*kbot == kwtop) {
940
941 /*        ==== 1-by-1 deflation window: not much to do ==== */
942
943         i__1 = kwtop;
944         i__2 = kwtop + kwtop * h_dim1;
945         sh[i__1].r = h__[i__2].r, sh[i__1].i = h__[i__2].i;
946         *ns = 1;
947         *nd = 0;
948 /* Computing MAX */
949         i__1 = kwtop + kwtop * h_dim1;
950         r__5 = smlnum, r__6 = ulp * ((r__1 = h__[i__1].r, abs(r__1)) + (r__2 =
951                  r_imag(&h__[kwtop + kwtop * h_dim1]), abs(r__2)));
952         if ((r__3 = s.r, abs(r__3)) + (r__4 = r_imag(&s), abs(r__4)) <= f2cmax(
953                 r__5,r__6)) {
954             *ns = 0;
955             *nd = 1;
956             if (kwtop > *ktop) {
957                 i__1 = kwtop + (kwtop - 1) * h_dim1;
958                 h__[i__1].r = 0.f, h__[i__1].i = 0.f;
959             }
960         }
961         work[1].r = 1.f, work[1].i = 0.f;
962         return 0;
963     }
964
965 /*     ==== Convert to spike-triangular form.  (In case of a */
966 /*     .    rare QR failure, this routine continues to do */
967 /*     .    aggressive early deflation using that part of */
968 /*     .    the deflation window that converged using INFQR */
969 /*     .    here and there to keep track.) ==== */
970
971     clacpy_("U", &jw, &jw, &h__[kwtop + kwtop * h_dim1], ldh, &t[t_offset], 
972             ldt);
973     i__1 = jw - 1;
974     i__2 = *ldh + 1;
975     i__3 = *ldt + 1;
976     ccopy_(&i__1, &h__[kwtop + 1 + kwtop * h_dim1], &i__2, &t[t_dim1 + 2], &
977             i__3);
978
979     claset_("A", &jw, &jw, &c_b1, &c_b2, &v[v_offset], ldv);
980     clahqr_(&c_true, &c_true, &jw, &c__1, &jw, &t[t_offset], ldt, &sh[kwtop], 
981             &c__1, &jw, &v[v_offset], ldv, &infqr);
982
983 /*     ==== Deflation detection loop ==== */
984
985     *ns = jw;
986     ilst = infqr + 1;
987     i__1 = jw;
988     for (knt = infqr + 1; knt <= i__1; ++knt) {
989
990 /*        ==== Small spike tip deflation test ==== */
991
992         i__2 = *ns + *ns * t_dim1;
993         foo = (r__1 = t[i__2].r, abs(r__1)) + (r__2 = r_imag(&t[*ns + *ns * 
994                 t_dim1]), abs(r__2));
995         if (foo == 0.f) {
996             foo = (r__1 = s.r, abs(r__1)) + (r__2 = r_imag(&s), abs(r__2));
997         }
998         i__2 = *ns * v_dim1 + 1;
999 /* Computing MAX */
1000         r__5 = smlnum, r__6 = ulp * foo;
1001         if (((r__1 = s.r, abs(r__1)) + (r__2 = r_imag(&s), abs(r__2))) * ((
1002                 r__3 = v[i__2].r, abs(r__3)) + (r__4 = r_imag(&v[*ns * v_dim1 
1003                 + 1]), abs(r__4))) <= f2cmax(r__5,r__6)) {
1004
1005 /*           ==== One more converged eigenvalue ==== */
1006
1007             --(*ns);
1008         } else {
1009
1010 /*           ==== One undeflatable eigenvalue.  Move it up out of the */
1011 /*           .    way.   (CTREXC can not fail in this case.) ==== */
1012
1013             ifst = *ns;
1014             ctrexc_("V", &jw, &t[t_offset], ldt, &v[v_offset], ldv, &ifst, &
1015                     ilst, &info);
1016             ++ilst;
1017         }
1018 /* L10: */
1019     }
1020
1021 /*        ==== Return to Hessenberg form ==== */
1022
1023     if (*ns == 0) {
1024         s.r = 0.f, s.i = 0.f;
1025     }
1026
1027     if (*ns < jw) {
1028
1029 /*        ==== sorting the diagonal of T improves accuracy for */
1030 /*        .    graded matrices.  ==== */
1031
1032         i__1 = *ns;
1033         for (i__ = infqr + 1; i__ <= i__1; ++i__) {
1034             ifst = i__;
1035             i__2 = *ns;
1036             for (j = i__ + 1; j <= i__2; ++j) {
1037                 i__3 = j + j * t_dim1;
1038                 i__4 = ifst + ifst * t_dim1;
1039                 if ((r__1 = t[i__3].r, abs(r__1)) + (r__2 = r_imag(&t[j + j * 
1040                         t_dim1]), abs(r__2)) > (r__3 = t[i__4].r, abs(r__3)) 
1041                         + (r__4 = r_imag(&t[ifst + ifst * t_dim1]), abs(r__4))
1042                         ) {
1043                     ifst = j;
1044                 }
1045 /* L20: */
1046             }
1047             ilst = i__;
1048             if (ifst != ilst) {
1049                 ctrexc_("V", &jw, &t[t_offset], ldt, &v[v_offset], ldv, &ifst,
1050                          &ilst, &info);
1051             }
1052 /* L30: */
1053         }
1054     }
1055
1056 /*     ==== Restore shift/eigenvalue array from T ==== */
1057
1058     i__1 = jw;
1059     for (i__ = infqr + 1; i__ <= i__1; ++i__) {
1060         i__2 = kwtop + i__ - 1;
1061         i__3 = i__ + i__ * t_dim1;
1062         sh[i__2].r = t[i__3].r, sh[i__2].i = t[i__3].i;
1063 /* L40: */
1064     }
1065
1066
1067     if (*ns < jw || s.r == 0.f && s.i == 0.f) {
1068         if (*ns > 1 && (s.r != 0.f || s.i != 0.f)) {
1069
1070 /*           ==== Reflect spike back into lower triangle ==== */
1071
1072             ccopy_(ns, &v[v_offset], ldv, &work[1], &c__1);
1073             i__1 = *ns;
1074             for (i__ = 1; i__ <= i__1; ++i__) {
1075                 i__2 = i__;
1076                 r_cnjg(&q__1, &work[i__]);
1077                 work[i__2].r = q__1.r, work[i__2].i = q__1.i;
1078 /* L50: */
1079             }
1080             beta.r = work[1].r, beta.i = work[1].i;
1081             clarfg_(ns, &beta, &work[2], &c__1, &tau);
1082             work[1].r = 1.f, work[1].i = 0.f;
1083
1084             i__1 = jw - 2;
1085             i__2 = jw - 2;
1086             claset_("L", &i__1, &i__2, &c_b1, &c_b1, &t[t_dim1 + 3], ldt);
1087
1088             r_cnjg(&q__1, &tau);
1089             clarf_("L", ns, &jw, &work[1], &c__1, &q__1, &t[t_offset], ldt, &
1090                     work[jw + 1]);
1091             clarf_("R", ns, ns, &work[1], &c__1, &tau, &t[t_offset], ldt, &
1092                     work[jw + 1]);
1093             clarf_("R", &jw, ns, &work[1], &c__1, &tau, &v[v_offset], ldv, &
1094                     work[jw + 1]);
1095
1096             i__1 = *lwork - jw;
1097             cgehrd_(&jw, &c__1, ns, &t[t_offset], ldt, &work[1], &work[jw + 1]
1098                     , &i__1, &info);
1099         }
1100
1101 /*        ==== Copy updated reduced window into place ==== */
1102
1103         if (kwtop > 1) {
1104             i__1 = kwtop + (kwtop - 1) * h_dim1;
1105             r_cnjg(&q__2, &v[v_dim1 + 1]);
1106             q__1.r = s.r * q__2.r - s.i * q__2.i, q__1.i = s.r * q__2.i + s.i 
1107                     * q__2.r;
1108             h__[i__1].r = q__1.r, h__[i__1].i = q__1.i;
1109         }
1110         clacpy_("U", &jw, &jw, &t[t_offset], ldt, &h__[kwtop + kwtop * h_dim1]
1111                 , ldh);
1112         i__1 = jw - 1;
1113         i__2 = *ldt + 1;
1114         i__3 = *ldh + 1;
1115         ccopy_(&i__1, &t[t_dim1 + 2], &i__2, &h__[kwtop + 1 + kwtop * h_dim1],
1116                  &i__3);
1117
1118 /*        ==== Accumulate orthogonal matrix in order update */
1119 /*        .    H and Z, if requested.  ==== */
1120
1121         if (*ns > 1 && (s.r != 0.f || s.i != 0.f)) {
1122             i__1 = *lwork - jw;
1123             cunmhr_("R", "N", &jw, ns, &c__1, ns, &t[t_offset], ldt, &work[1],
1124                      &v[v_offset], ldv, &work[jw + 1], &i__1, &info);
1125         }
1126
1127 /*        ==== Update vertical slab in H ==== */
1128
1129         if (*wantt) {
1130             ltop = 1;
1131         } else {
1132             ltop = *ktop;
1133         }
1134         i__1 = kwtop - 1;
1135         i__2 = *nv;
1136         for (krow = ltop; i__2 < 0 ? krow >= i__1 : krow <= i__1; krow += 
1137                 i__2) {
1138 /* Computing MIN */
1139             i__3 = *nv, i__4 = kwtop - krow;
1140             kln = f2cmin(i__3,i__4);
1141             cgemm_("N", "N", &kln, &jw, &jw, &c_b2, &h__[krow + kwtop * 
1142                     h_dim1], ldh, &v[v_offset], ldv, &c_b1, &wv[wv_offset], 
1143                     ldwv);
1144             clacpy_("A", &kln, &jw, &wv[wv_offset], ldwv, &h__[krow + kwtop * 
1145                     h_dim1], ldh);
1146 /* L60: */
1147         }
1148
1149 /*        ==== Update horizontal slab in H ==== */
1150
1151         if (*wantt) {
1152             i__2 = *n;
1153             i__1 = *nh;
1154             for (kcol = *kbot + 1; i__1 < 0 ? kcol >= i__2 : kcol <= i__2; 
1155                     kcol += i__1) {
1156 /* Computing MIN */
1157                 i__3 = *nh, i__4 = *n - kcol + 1;
1158                 kln = f2cmin(i__3,i__4);
1159                 cgemm_("C", "N", &jw, &kln, &jw, &c_b2, &v[v_offset], ldv, &
1160                         h__[kwtop + kcol * h_dim1], ldh, &c_b1, &t[t_offset], 
1161                         ldt);
1162                 clacpy_("A", &jw, &kln, &t[t_offset], ldt, &h__[kwtop + kcol *
1163                          h_dim1], ldh);
1164 /* L70: */
1165             }
1166         }
1167
1168 /*        ==== Update vertical slab in Z ==== */
1169
1170         if (*wantz) {
1171             i__1 = *ihiz;
1172             i__2 = *nv;
1173             for (krow = *iloz; i__2 < 0 ? krow >= i__1 : krow <= i__1; krow +=
1174                      i__2) {
1175 /* Computing MIN */
1176                 i__3 = *nv, i__4 = *ihiz - krow + 1;
1177                 kln = f2cmin(i__3,i__4);
1178                 cgemm_("N", "N", &kln, &jw, &jw, &c_b2, &z__[krow + kwtop * 
1179                         z_dim1], ldz, &v[v_offset], ldv, &c_b1, &wv[wv_offset]
1180                         , ldwv);
1181                 clacpy_("A", &kln, &jw, &wv[wv_offset], ldwv, &z__[krow + 
1182                         kwtop * z_dim1], ldz);
1183 /* L80: */
1184             }
1185         }
1186     }
1187
1188 /*     ==== Return the number of deflations ... ==== */
1189
1190     *nd = jw - *ns;
1191
1192 /*     ==== ... and the number of shifts. (Subtracting */
1193 /*     .    INFQR from the spike length takes care */
1194 /*     .    of the case of a rare QR failure while */
1195 /*     .    calculating eigenvalues of the deflation */
1196 /*     .    window.)  ==== */
1197
1198     *ns -= infqr;
1199
1200 /*      ==== Return optimal workspace. ==== */
1201
1202     r__1 = (real) lwkopt;
1203     q__1.r = r__1, q__1.i = 0.f;
1204     work[1].r = q__1.r, work[1].i = q__1.i;
1205
1206 /*     ==== End of CLAQR2 ==== */
1207
1208     return 0;
1209 } /* claqr2_ */
1210