C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / clahr2.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
519 /* > \brief \b CLAHR2 reduces the specified number of first columns of a general rectangular matrix A so that 
520 elements below the specified subdiagonal are zero, and returns auxiliary matrices which are needed to 
521 apply the transformation to the unreduced part */
522 /* of A. */
523
524 /*  =========== DOCUMENTATION =========== */
525
526 /* Online html documentation available at */
527 /*            http://www.netlib.org/lapack/explore-html/ */
528
529 /* > \htmlonly */
530 /* > Download CLAHR2 + dependencies */
531 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/clahr2.
532 f"> */
533 /* > [TGZ]</a> */
534 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/clahr2.
535 f"> */
536 /* > [ZIP]</a> */
537 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/clahr2.
538 f"> */
539 /* > [TXT]</a> */
540 /* > \endhtmlonly */
541
542 /*  Definition: */
543 /*  =========== */
544
545 /*       SUBROUTINE CLAHR2( N, K, NB, A, LDA, TAU, T, LDT, Y, LDY ) */
546
547 /*       INTEGER            K, LDA, LDT, LDY, N, NB */
548 /*       COMPLEX            A( LDA, * ), T( LDT, NB ), TAU( NB ), */
549 /*      $                   Y( LDY, NB ) */
550
551
552 /* > \par Purpose: */
553 /*  ============= */
554 /* > */
555 /* > \verbatim */
556 /* > */
557 /* > CLAHR2 reduces the first NB columns of A complex general n-BY-(n-k+1) */
558 /* > matrix A so that elements below the k-th subdiagonal are zero. The */
559 /* > reduction is performed by an unitary similarity transformation */
560 /* > Q**H * A * Q. The routine returns the matrices V and T which determine */
561 /* > Q as a block reflector I - V*T*v**H, and also the matrix Y = A * V * T. */
562 /* > */
563 /* > This is an auxiliary routine called by CGEHRD. */
564 /* > \endverbatim */
565
566 /*  Arguments: */
567 /*  ========== */
568
569 /* > \param[in] N */
570 /* > \verbatim */
571 /* >          N is INTEGER */
572 /* >          The order of the matrix A. */
573 /* > \endverbatim */
574 /* > */
575 /* > \param[in] K */
576 /* > \verbatim */
577 /* >          K is INTEGER */
578 /* >          The offset for the reduction. Elements below the k-th */
579 /* >          subdiagonal in the first NB columns are reduced to zero. */
580 /* >          K < N. */
581 /* > \endverbatim */
582 /* > */
583 /* > \param[in] NB */
584 /* > \verbatim */
585 /* >          NB is INTEGER */
586 /* >          The number of columns to be reduced. */
587 /* > \endverbatim */
588 /* > */
589 /* > \param[in,out] A */
590 /* > \verbatim */
591 /* >          A is COMPLEX array, dimension (LDA,N-K+1) */
592 /* >          On entry, the n-by-(n-k+1) general matrix A. */
593 /* >          On exit, the elements on and above the k-th subdiagonal in */
594 /* >          the first NB columns are overwritten with the corresponding */
595 /* >          elements of the reduced matrix; the elements below the k-th */
596 /* >          subdiagonal, with the array TAU, represent the matrix Q as a */
597 /* >          product of elementary reflectors. The other columns of A are */
598 /* >          unchanged. See Further Details. */
599 /* > \endverbatim */
600 /* > */
601 /* > \param[in] LDA */
602 /* > \verbatim */
603 /* >          LDA is INTEGER */
604 /* >          The leading dimension of the array A.  LDA >= f2cmax(1,N). */
605 /* > \endverbatim */
606 /* > */
607 /* > \param[out] TAU */
608 /* > \verbatim */
609 /* >          TAU is COMPLEX array, dimension (NB) */
610 /* >          The scalar factors of the elementary reflectors. See Further */
611 /* >          Details. */
612 /* > \endverbatim */
613 /* > */
614 /* > \param[out] T */
615 /* > \verbatim */
616 /* >          T is COMPLEX array, dimension (LDT,NB) */
617 /* >          The upper triangular matrix T. */
618 /* > \endverbatim */
619 /* > */
620 /* > \param[in] LDT */
621 /* > \verbatim */
622 /* >          LDT is INTEGER */
623 /* >          The leading dimension of the array T.  LDT >= NB. */
624 /* > \endverbatim */
625 /* > */
626 /* > \param[out] Y */
627 /* > \verbatim */
628 /* >          Y is COMPLEX array, dimension (LDY,NB) */
629 /* >          The n-by-nb matrix Y. */
630 /* > \endverbatim */
631 /* > */
632 /* > \param[in] LDY */
633 /* > \verbatim */
634 /* >          LDY is INTEGER */
635 /* >          The leading dimension of the array Y. LDY >= N. */
636 /* > \endverbatim */
637
638 /*  Authors: */
639 /*  ======== */
640
641 /* > \author Univ. of Tennessee */
642 /* > \author Univ. of California Berkeley */
643 /* > \author Univ. of Colorado Denver */
644 /* > \author NAG Ltd. */
645
646 /* > \date December 2016 */
647
648 /* > \ingroup complexOTHERauxiliary */
649
650 /* > \par Further Details: */
651 /*  ===================== */
652 /* > */
653 /* > \verbatim */
654 /* > */
655 /* >  The matrix Q is represented as a product of nb elementary reflectors */
656 /* > */
657 /* >     Q = H(1) H(2) . . . H(nb). */
658 /* > */
659 /* >  Each H(i) has the form */
660 /* > */
661 /* >     H(i) = I - tau * v * v**H */
662 /* > */
663 /* >  where tau is a complex scalar, and v is a complex vector with */
664 /* >  v(1:i+k-1) = 0, v(i+k) = 1; v(i+k+1:n) is stored on exit in */
665 /* >  A(i+k+1:n,i), and tau in TAU(i). */
666 /* > */
667 /* >  The elements of the vectors v together form the (n-k+1)-by-nb matrix */
668 /* >  V which is needed, with T and Y, to apply the transformation to the */
669 /* >  unreduced part of the matrix, using an update of the form: */
670 /* >  A := (I - V*T*V**H) * (A - Y*V**H). */
671 /* > */
672 /* >  The contents of A on exit are illustrated by the following example */
673 /* >  with n = 7, k = 3 and nb = 2: */
674 /* > */
675 /* >     ( a   a   a   a   a ) */
676 /* >     ( a   a   a   a   a ) */
677 /* >     ( a   a   a   a   a ) */
678 /* >     ( h   h   a   a   a ) */
679 /* >     ( v1  h   a   a   a ) */
680 /* >     ( v1  v2  a   a   a ) */
681 /* >     ( v1  v2  a   a   a ) */
682 /* > */
683 /* >  where a denotes an element of the original matrix A, h denotes a */
684 /* >  modified element of the upper Hessenberg matrix H, and vi denotes an */
685 /* >  element of the vector defining H(i). */
686 /* > */
687 /* >  This subroutine is a slight modification of LAPACK-3.0's DLAHRD */
688 /* >  incorporating improvements proposed by Quintana-Orti and Van de */
689 /* >  Gejin. Note that the entries of A(1:K,2:NB) differ from those */
690 /* >  returned by the original LAPACK-3.0's DLAHRD routine. (This */
691 /* >  subroutine is not backward compatible with LAPACK-3.0's DLAHRD.) */
692 /* > \endverbatim */
693
694 /* > \par References: */
695 /*  ================ */
696 /* > */
697 /* >  Gregorio Quintana-Orti and Robert van de Geijn, "Improving the */
698 /* >  performance of reduction to Hessenberg form," ACM Transactions on */
699 /* >  Mathematical Software, 32(2):180-194, June 2006. */
700 /* > */
701 /*  ===================================================================== */
702 /* Subroutine */ int clahr2_(integer *n, integer *k, integer *nb, complex *a, 
703         integer *lda, complex *tau, complex *t, integer *ldt, complex *y, 
704         integer *ldy)
705 {
706     /* System generated locals */
707     integer a_dim1, a_offset, t_dim1, t_offset, y_dim1, y_offset, i__1, i__2, 
708             i__3;
709     complex q__1;
710
711     /* Local variables */
712     integer i__;
713     extern /* Subroutine */ int cscal_(integer *, complex *, complex *, 
714             integer *), cgemm_(char *, char *, integer *, integer *, integer *
715             , complex *, complex *, integer *, complex *, integer *, complex *
716             , complex *, integer *), cgemv_(char *, integer *,
717              integer *, complex *, complex *, integer *, complex *, integer *,
718              complex *, complex *, integer *), ccopy_(integer *, 
719             complex *, integer *, complex *, integer *), ctrmm_(char *, char *
720             , char *, char *, integer *, integer *, complex *, complex *, 
721             integer *, complex *, integer *), 
722             caxpy_(integer *, complex *, complex *, integer *, complex *, 
723             integer *), ctrmv_(char *, char *, char *, integer *, complex *, 
724             integer *, complex *, integer *);
725     complex ei;
726     extern /* Subroutine */ int clarfg_(integer *, complex *, complex *, 
727             integer *, complex *), clacgv_(integer *, complex *, integer *), 
728             clacpy_(char *, integer *, integer *, complex *, integer *, 
729             complex *, integer *);
730
731
732 /*  -- LAPACK auxiliary routine (version 3.7.0) -- */
733 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
734 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
735 /*     December 2016 */
736
737
738 /*  ===================================================================== */
739
740
741 /*     Quick return if possible */
742
743     /* Parameter adjustments */
744     --tau;
745     a_dim1 = *lda;
746     a_offset = 1 + a_dim1 * 1;
747     a -= a_offset;
748     t_dim1 = *ldt;
749     t_offset = 1 + t_dim1 * 1;
750     t -= t_offset;
751     y_dim1 = *ldy;
752     y_offset = 1 + y_dim1 * 1;
753     y -= y_offset;
754
755     /* Function Body */
756     if (*n <= 1) {
757         return 0;
758     }
759
760     i__1 = *nb;
761     for (i__ = 1; i__ <= i__1; ++i__) {
762         if (i__ > 1) {
763
764 /*           Update A(K+1:N,I) */
765
766 /*           Update I-th column of A - Y * V**H */
767
768             i__2 = i__ - 1;
769             clacgv_(&i__2, &a[*k + i__ - 1 + a_dim1], lda);
770             i__2 = *n - *k;
771             i__3 = i__ - 1;
772             q__1.r = -1.f, q__1.i = 0.f;
773             cgemv_("NO TRANSPOSE", &i__2, &i__3, &q__1, &y[*k + 1 + y_dim1], 
774                     ldy, &a[*k + i__ - 1 + a_dim1], lda, &c_b2, &a[*k + 1 + 
775                     i__ * a_dim1], &c__1);
776             i__2 = i__ - 1;
777             clacgv_(&i__2, &a[*k + i__ - 1 + a_dim1], lda);
778
779 /*           Apply I - V * T**H * V**H to this column (call it b) from the */
780 /*           left, using the last column of T as workspace */
781
782 /*           Let  V = ( V1 )   and   b = ( b1 )   (first I-1 rows) */
783 /*                    ( V2 )             ( b2 ) */
784
785 /*           where V1 is unit lower triangular */
786
787 /*           w := V1**H * b1 */
788
789             i__2 = i__ - 1;
790             ccopy_(&i__2, &a[*k + 1 + i__ * a_dim1], &c__1, &t[*nb * t_dim1 + 
791                     1], &c__1);
792             i__2 = i__ - 1;
793             ctrmv_("Lower", "Conjugate transpose", "UNIT", &i__2, &a[*k + 1 + 
794                     a_dim1], lda, &t[*nb * t_dim1 + 1], &c__1);
795
796 /*           w := w + V2**H * b2 */
797
798             i__2 = *n - *k - i__ + 1;
799             i__3 = i__ - 1;
800             cgemv_("Conjugate transpose", &i__2, &i__3, &c_b2, &a[*k + i__ + 
801                     a_dim1], lda, &a[*k + i__ + i__ * a_dim1], &c__1, &c_b2, &
802                     t[*nb * t_dim1 + 1], &c__1);
803
804 /*           w := T**H * w */
805
806             i__2 = i__ - 1;
807             ctrmv_("Upper", "Conjugate transpose", "NON-UNIT", &i__2, &t[
808                     t_offset], ldt, &t[*nb * t_dim1 + 1], &c__1);
809
810 /*           b2 := b2 - V2*w */
811
812             i__2 = *n - *k - i__ + 1;
813             i__3 = i__ - 1;
814             q__1.r = -1.f, q__1.i = 0.f;
815             cgemv_("NO TRANSPOSE", &i__2, &i__3, &q__1, &a[*k + i__ + a_dim1],
816                      lda, &t[*nb * t_dim1 + 1], &c__1, &c_b2, &a[*k + i__ + 
817                     i__ * a_dim1], &c__1);
818
819 /*           b1 := b1 - V1*w */
820
821             i__2 = i__ - 1;
822             ctrmv_("Lower", "NO TRANSPOSE", "UNIT", &i__2, &a[*k + 1 + a_dim1]
823                     , lda, &t[*nb * t_dim1 + 1], &c__1);
824             i__2 = i__ - 1;
825             q__1.r = -1.f, q__1.i = 0.f;
826             caxpy_(&i__2, &q__1, &t[*nb * t_dim1 + 1], &c__1, &a[*k + 1 + i__ 
827                     * a_dim1], &c__1);
828
829             i__2 = *k + i__ - 1 + (i__ - 1) * a_dim1;
830             a[i__2].r = ei.r, a[i__2].i = ei.i;
831         }
832
833 /*        Generate the elementary reflector H(I) to annihilate */
834 /*        A(K+I+1:N,I) */
835
836         i__2 = *n - *k - i__ + 1;
837 /* Computing MIN */
838         i__3 = *k + i__ + 1;
839         clarfg_(&i__2, &a[*k + i__ + i__ * a_dim1], &a[f2cmin(i__3,*n) + i__ * 
840                 a_dim1], &c__1, &tau[i__]);
841         i__2 = *k + i__ + i__ * a_dim1;
842         ei.r = a[i__2].r, ei.i = a[i__2].i;
843         i__2 = *k + i__ + i__ * a_dim1;
844         a[i__2].r = 1.f, a[i__2].i = 0.f;
845
846 /*        Compute  Y(K+1:N,I) */
847
848         i__2 = *n - *k;
849         i__3 = *n - *k - i__ + 1;
850         cgemv_("NO TRANSPOSE", &i__2, &i__3, &c_b2, &a[*k + 1 + (i__ + 1) * 
851                 a_dim1], lda, &a[*k + i__ + i__ * a_dim1], &c__1, &c_b1, &y[*
852                 k + 1 + i__ * y_dim1], &c__1);
853         i__2 = *n - *k - i__ + 1;
854         i__3 = i__ - 1;
855         cgemv_("Conjugate transpose", &i__2, &i__3, &c_b2, &a[*k + i__ + 
856                 a_dim1], lda, &a[*k + i__ + i__ * a_dim1], &c__1, &c_b1, &t[
857                 i__ * t_dim1 + 1], &c__1);
858         i__2 = *n - *k;
859         i__3 = i__ - 1;
860         q__1.r = -1.f, q__1.i = 0.f;
861         cgemv_("NO TRANSPOSE", &i__2, &i__3, &q__1, &y[*k + 1 + y_dim1], ldy, 
862                 &t[i__ * t_dim1 + 1], &c__1, &c_b2, &y[*k + 1 + i__ * y_dim1],
863                  &c__1);
864         i__2 = *n - *k;
865         cscal_(&i__2, &tau[i__], &y[*k + 1 + i__ * y_dim1], &c__1);
866
867 /*        Compute T(1:I,I) */
868
869         i__2 = i__ - 1;
870         i__3 = i__;
871         q__1.r = -tau[i__3].r, q__1.i = -tau[i__3].i;
872         cscal_(&i__2, &q__1, &t[i__ * t_dim1 + 1], &c__1);
873         i__2 = i__ - 1;
874         ctrmv_("Upper", "No Transpose", "NON-UNIT", &i__2, &t[t_offset], ldt, 
875                 &t[i__ * t_dim1 + 1], &c__1)
876                 ;
877         i__2 = i__ + i__ * t_dim1;
878         i__3 = i__;
879         t[i__2].r = tau[i__3].r, t[i__2].i = tau[i__3].i;
880
881 /* L10: */
882     }
883     i__1 = *k + *nb + *nb * a_dim1;
884     a[i__1].r = ei.r, a[i__1].i = ei.i;
885
886 /*     Compute Y(1:K,1:NB) */
887
888     clacpy_("ALL", k, nb, &a[(a_dim1 << 1) + 1], lda, &y[y_offset], ldy);
889     ctrmm_("RIGHT", "Lower", "NO TRANSPOSE", "UNIT", k, nb, &c_b2, &a[*k + 1 
890             + a_dim1], lda, &y[y_offset], ldy);
891     if (*n > *k + *nb) {
892         i__1 = *n - *k - *nb;
893         cgemm_("NO TRANSPOSE", "NO TRANSPOSE", k, nb, &i__1, &c_b2, &a[(*nb + 
894                 2) * a_dim1 + 1], lda, &a[*k + 1 + *nb + a_dim1], lda, &c_b2, 
895                 &y[y_offset], ldy);
896     }
897     ctrmm_("RIGHT", "Upper", "NO TRANSPOSE", "NON-UNIT", k, nb, &c_b2, &t[
898             t_offset], ldt, &y[y_offset], ldy);
899
900     return 0;
901
902 /*     End of CLAHR2 */
903
904 } /* clahr2_ */
905