C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / dlasyf_aa.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 doublereal c_b6 = -1.;
516 static integer c__1 = 1;
517 static doublereal c_b8 = 1.;
518 static doublereal c_b22 = 0.;
519
520 /* > \brief \b DLASYF_AA */
521
522 /*  =========== DOCUMENTATION =========== */
523
524 /* Online html documentation available at */
525 /*            http://www.netlib.org/lapack/explore-html/ */
526
527 /* > \htmlonly */
528 /* > Download DLASYF_AA + dependencies */
529 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlasyf_
530 aa.f"> */
531 /* > [TGZ]</a> */
532 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlasyf_
533 aa.f"> */
534 /* > [ZIP]</a> */
535 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlasyf_
536 aa.f"> */
537 /* > [TXT]</a> */
538 /* > \endhtmlonly */
539
540 /*  Definition: */
541 /*  =========== */
542
543 /*       SUBROUTINE DLASYF_AA( UPLO, J1, M, NB, A, LDA, IPIV, */
544 /*                             H, LDH, WORK ) */
545
546 /*       CHARACTER          UPLO */
547 /*       INTEGER            J1, M, NB, LDA, LDH */
548 /*       INTEGER            IPIV( * ) */
549 /*       DOUBLE PRECISION   A( LDA, * ), H( LDH, * ), WORK( * ) */
550
551
552 /* > \par Purpose: */
553 /*  ============= */
554 /* > */
555 /* > \verbatim */
556 /* > */
557 /* > DLATRF_AA factorizes a panel of a real symmetric matrix A using */
558 /* > the Aasen's algorithm. The panel consists of a set of NB rows of A */
559 /* > when UPLO is U, or a set of NB columns when UPLO is L. */
560 /* > */
561 /* > In order to factorize the panel, the Aasen's algorithm requires the */
562 /* > last row, or column, of the previous panel. The first row, or column, */
563 /* > of A is set to be the first row, or column, of an identity matrix, */
564 /* > which is used to factorize the first panel. */
565 /* > */
566 /* > The resulting J-th row of U, or J-th column of L, is stored in the */
567 /* > (J-1)-th row, or column, of A (without the unit diagonals), while */
568 /* > the diagonal and subdiagonal of A are overwritten by those of T. */
569 /* > */
570 /* > \endverbatim */
571
572 /*  Arguments: */
573 /*  ========== */
574
575 /* > \param[in] UPLO */
576 /* > \verbatim */
577 /* >          UPLO is CHARACTER*1 */
578 /* >          = 'U':  Upper triangle of A is stored; */
579 /* >          = 'L':  Lower triangle of A is stored. */
580 /* > \endverbatim */
581 /* > */
582 /* > \param[in] J1 */
583 /* > \verbatim */
584 /* >          J1 is INTEGER */
585 /* >          The location of the first row, or column, of the panel */
586 /* >          within the submatrix of A, passed to this routine, e.g., */
587 /* >          when called by DSYTRF_AA, for the first panel, J1 is 1, */
588 /* >          while for the remaining panels, J1 is 2. */
589 /* > \endverbatim */
590 /* > */
591 /* > \param[in] M */
592 /* > \verbatim */
593 /* >          M is INTEGER */
594 /* >          The dimension of the submatrix. M >= 0. */
595 /* > \endverbatim */
596 /* > */
597 /* > \param[in] NB */
598 /* > \verbatim */
599 /* >          NB is INTEGER */
600 /* >          The dimension of the panel to be facotorized. */
601 /* > \endverbatim */
602 /* > */
603 /* > \param[in,out] A */
604 /* > \verbatim */
605 /* >          A is DOUBLE PRECISION array, dimension (LDA,M) for */
606 /* >          the first panel, while dimension (LDA,M+1) for the */
607 /* >          remaining panels. */
608 /* > */
609 /* >          On entry, A contains the last row, or column, of */
610 /* >          the previous panel, and the trailing submatrix of A */
611 /* >          to be factorized, except for the first panel, only */
612 /* >          the panel is passed. */
613 /* > */
614 /* >          On exit, the leading panel is factorized. */
615 /* > \endverbatim */
616 /* > */
617 /* > \param[in] LDA */
618 /* > \verbatim */
619 /* >          LDA is INTEGER */
620 /* >          The leading dimension of the array A.  LDA >= f2cmax(1,M). */
621 /* > \endverbatim */
622 /* > */
623 /* > \param[out] IPIV */
624 /* > \verbatim */
625 /* >          IPIV is INTEGER array, dimension (M) */
626 /* >          Details of the row and column interchanges, */
627 /* >          the row and column k were interchanged with the row and */
628 /* >          column IPIV(k). */
629 /* > \endverbatim */
630 /* > */
631 /* > \param[in,out] H */
632 /* > \verbatim */
633 /* >          H is DOUBLE PRECISION workspace, dimension (LDH,NB). */
634 /* > */
635 /* > \endverbatim */
636 /* > */
637 /* > \param[in] LDH */
638 /* > \verbatim */
639 /* >          LDH is INTEGER */
640 /* >          The leading dimension of the workspace H. LDH >= f2cmax(1,M). */
641 /* > \endverbatim */
642 /* > */
643 /* > \param[out] WORK */
644 /* > \verbatim */
645 /* >          WORK is DOUBLE PRECISION workspace, dimension (M). */
646 /* > \endverbatim */
647 /* > */
648
649 /*  Authors: */
650 /*  ======== */
651
652 /* > \author Univ. of Tennessee */
653 /* > \author Univ. of California Berkeley */
654 /* > \author Univ. of Colorado Denver */
655 /* > \author NAG Ltd. */
656
657 /* > \date November 2017 */
658
659 /* > \ingroup doubleSYcomputational */
660
661 /*  ===================================================================== */
662 /* Subroutine */ int dlasyf_aa_(char *uplo, integer *j1, integer *m, integer 
663         *nb, doublereal *a, integer *lda, integer *ipiv, doublereal *h__, 
664         integer *ldh, doublereal *work)
665 {
666     /* System generated locals */
667     integer a_dim1, a_offset, h_dim1, h_offset, i__1;
668
669     /* Local variables */
670     integer j, k;
671     doublereal alpha;
672     extern /* Subroutine */ int dscal_(integer *, doublereal *, doublereal *, 
673             integer *);
674     extern logical lsame_(char *, char *);
675     extern /* Subroutine */ int dgemv_(char *, integer *, integer *, 
676             doublereal *, doublereal *, integer *, doublereal *, integer *, 
677             doublereal *, doublereal *, integer *), dcopy_(integer *, 
678             doublereal *, integer *, doublereal *, integer *), dswap_(integer 
679             *, doublereal *, integer *, doublereal *, integer *), daxpy_(
680             integer *, doublereal *, doublereal *, integer *, doublereal *, 
681             integer *);
682     integer i1, k1, i2, mj;
683     extern integer idamax_(integer *, doublereal *, integer *);
684     extern /* Subroutine */ int dlaset_(char *, integer *, integer *, 
685             doublereal *, doublereal *, doublereal *, integer *);
686     doublereal piv;
687
688
689 /*  -- LAPACK computational routine (version 3.8.0) -- */
690 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
691 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
692 /*     November 2017 */
693
694
695
696 /*  ===================================================================== */
697
698
699     /* Parameter adjustments */
700     a_dim1 = *lda;
701     a_offset = 1 + a_dim1 * 1;
702     a -= a_offset;
703     --ipiv;
704     h_dim1 = *ldh;
705     h_offset = 1 + h_dim1 * 1;
706     h__ -= h_offset;
707     --work;
708
709     /* Function Body */
710     j = 1;
711
712 /*     K1 is the first column of the panel to be factorized */
713 /*     i.e.,  K1 is 2 for the first block column, and 1 for the rest of the blocks */
714
715     k1 = 2 - *j1 + 1;
716
717     if (lsame_(uplo, "U")) {
718
719 /*        ..................................................... */
720 /*        Factorize A as U**T*D*U using the upper triangle of A */
721 /*        ..................................................... */
722
723 L10:
724         if (j > f2cmin(*m,*nb)) {
725             goto L20;
726         }
727
728 /*        K is the column to be factorized */
729 /*         when being called from DSYTRF_AA, */
730 /*         > for the first block column, J1 is 1, hence J1+J-1 is J, */
731 /*         > for the rest of the columns, J1 is 2, and J1+J-1 is J+1, */
732
733         k = *j1 + j - 1;
734         if (j == *m) {
735
736 /*            Only need to compute T(J, J) */
737
738             mj = 1;
739         } else {
740             mj = *m - j + 1;
741         }
742
743 /*        H(J:M, J) := A(J, J:M) - H(J:M, 1:(J-1)) * L(J1:(J-1), J), */
744 /*         where H(J:M, J) has been initialized to be A(J, J:M) */
745
746         if (k > 2) {
747
748 /*        K is the column to be factorized */
749 /*         > for the first block column, K is J, skipping the first two */
750 /*           columns */
751 /*         > for the rest of the columns, K is J+1, skipping only the */
752 /*           first column */
753
754             i__1 = j - k1;
755             dgemv_("No transpose", &mj, &i__1, &c_b6, &h__[j + k1 * h_dim1], 
756                     ldh, &a[j * a_dim1 + 1], &c__1, &c_b8, &h__[j + j * 
757                     h_dim1], &c__1);
758         }
759
760 /*        Copy H(i:M, i) into WORK */
761
762         dcopy_(&mj, &h__[j + j * h_dim1], &c__1, &work[1], &c__1);
763
764         if (j > k1) {
765
766 /*           Compute WORK := WORK - L(J-1, J:M) * T(J-1,J), */
767 /*            where A(J-1, J) stores T(J-1, J) and A(J-2, J:M) stores U(J-1, J:M) */
768
769             alpha = -a[k - 1 + j * a_dim1];
770             daxpy_(&mj, &alpha, &a[k - 2 + j * a_dim1], lda, &work[1], &c__1);
771         }
772
773 /*        Set A(J, J) = T(J, J) */
774
775         a[k + j * a_dim1] = work[1];
776
777         if (j < *m) {
778
779 /*           Compute WORK(2:M) = T(J, J) L(J, (J+1):M) */
780 /*            where A(J, J) stores T(J, J) and A(J-1, (J+1):M) stores U(J, (J+1):M) */
781
782             if (k > 1) {
783                 alpha = -a[k + j * a_dim1];
784                 i__1 = *m - j;
785                 daxpy_(&i__1, &alpha, &a[k - 1 + (j + 1) * a_dim1], lda, &
786                         work[2], &c__1);
787             }
788
789 /*           Find f2cmax(|WORK(2:M)|) */
790
791             i__1 = *m - j;
792             i2 = idamax_(&i__1, &work[2], &c__1) + 1;
793             piv = work[i2];
794
795 /*           Apply symmetric pivot */
796
797             if (i2 != 2 && piv != 0.) {
798
799 /*              Swap WORK(I1) and WORK(I2) */
800
801                 i1 = 2;
802                 work[i2] = work[i1];
803                 work[i1] = piv;
804
805 /*              Swap A(I1, I1+1:M) with A(I1+1:M, I2) */
806
807                 i1 = i1 + j - 1;
808                 i2 = i2 + j - 1;
809                 i__1 = i2 - i1 - 1;
810                 dswap_(&i__1, &a[*j1 + i1 - 1 + (i1 + 1) * a_dim1], lda, &a[*
811                         j1 + i1 + i2 * a_dim1], &c__1);
812
813 /*              Swap A(I1, I2+1:M) with A(I2, I2+1:M) */
814
815                 if (i2 < *m) {
816                     i__1 = *m - i2;
817                     dswap_(&i__1, &a[*j1 + i1 - 1 + (i2 + 1) * a_dim1], lda, &
818                             a[*j1 + i2 - 1 + (i2 + 1) * a_dim1], lda);
819                 }
820
821 /*              Swap A(I1, I1) with A(I2,I2) */
822
823                 piv = a[i1 + *j1 - 1 + i1 * a_dim1];
824                 a[*j1 + i1 - 1 + i1 * a_dim1] = a[*j1 + i2 - 1 + i2 * a_dim1];
825                 a[*j1 + i2 - 1 + i2 * a_dim1] = piv;
826
827 /*              Swap H(I1, 1:J1) with H(I2, 1:J1) */
828
829                 i__1 = i1 - 1;
830                 dswap_(&i__1, &h__[i1 + h_dim1], ldh, &h__[i2 + h_dim1], ldh);
831                 ipiv[i1] = i2;
832
833                 if (i1 > k1 - 1) {
834
835 /*                 Swap L(1:I1-1, I1) with L(1:I1-1, I2), */
836 /*                  skipping the first column */
837
838                     i__1 = i1 - k1 + 1;
839                     dswap_(&i__1, &a[i1 * a_dim1 + 1], &c__1, &a[i2 * a_dim1 
840                             + 1], &c__1);
841                 }
842             } else {
843                 ipiv[j + 1] = j + 1;
844             }
845
846 /*           Set A(J, J+1) = T(J, J+1) */
847
848             a[k + (j + 1) * a_dim1] = work[2];
849
850             if (j < *nb) {
851
852 /*              Copy A(J+1:M, J+1) into H(J:M, J), */
853
854                 i__1 = *m - j;
855                 dcopy_(&i__1, &a[k + 1 + (j + 1) * a_dim1], lda, &h__[j + 1 + 
856                         (j + 1) * h_dim1], &c__1);
857             }
858
859 /*           Compute L(J+2, J+1) = WORK( 3:M ) / T(J, J+1), */
860 /*            where A(J, J+1) = T(J, J+1) and A(J+2:M, J) = L(J+2:M, J+1) */
861
862             if (j < *m - 1) {
863                 if (a[k + (j + 1) * a_dim1] != 0.) {
864                     alpha = 1. / a[k + (j + 1) * a_dim1];
865                     i__1 = *m - j - 1;
866                     dcopy_(&i__1, &work[3], &c__1, &a[k + (j + 2) * a_dim1], 
867                             lda);
868                     i__1 = *m - j - 1;
869                     dscal_(&i__1, &alpha, &a[k + (j + 2) * a_dim1], lda);
870                 } else {
871                     i__1 = *m - j - 1;
872                     dlaset_("Full", &c__1, &i__1, &c_b22, &c_b22, &a[k + (j + 
873                             2) * a_dim1], lda);
874                 }
875             }
876         }
877         ++j;
878         goto L10;
879 L20:
880
881         ;
882     } else {
883
884 /*        ..................................................... */
885 /*        Factorize A as L*D*L**T using the lower triangle of A */
886 /*        ..................................................... */
887
888 L30:
889         if (j > f2cmin(*m,*nb)) {
890             goto L40;
891         }
892
893 /*        K is the column to be factorized */
894 /*         when being called from DSYTRF_AA, */
895 /*         > for the first block column, J1 is 1, hence J1+J-1 is J, */
896 /*         > for the rest of the columns, J1 is 2, and J1+J-1 is J+1, */
897
898         k = *j1 + j - 1;
899         if (j == *m) {
900
901 /*            Only need to compute T(J, J) */
902
903             mj = 1;
904         } else {
905             mj = *m - j + 1;
906         }
907
908 /*        H(J:M, J) := A(J:M, J) - H(J:M, 1:(J-1)) * L(J, J1:(J-1))^T, */
909 /*         where H(J:M, J) has been initialized to be A(J:M, J) */
910
911         if (k > 2) {
912
913 /*        K is the column to be factorized */
914 /*         > for the first block column, K is J, skipping the first two */
915 /*           columns */
916 /*         > for the rest of the columns, K is J+1, skipping only the */
917 /*           first column */
918
919             i__1 = j - k1;
920             dgemv_("No transpose", &mj, &i__1, &c_b6, &h__[j + k1 * h_dim1], 
921                     ldh, &a[j + a_dim1], lda, &c_b8, &h__[j + j * h_dim1], &
922                     c__1);
923         }
924
925 /*        Copy H(J:M, J) into WORK */
926
927         dcopy_(&mj, &h__[j + j * h_dim1], &c__1, &work[1], &c__1);
928
929         if (j > k1) {
930
931 /*           Compute WORK := WORK - L(J:M, J-1) * T(J-1,J), */
932 /*            where A(J-1, J) = T(J-1, J) and A(J, J-2) = L(J, J-1) */
933
934             alpha = -a[j + (k - 1) * a_dim1];
935             daxpy_(&mj, &alpha, &a[j + (k - 2) * a_dim1], &c__1, &work[1], &
936                     c__1);
937         }
938
939 /*        Set A(J, J) = T(J, J) */
940
941         a[j + k * a_dim1] = work[1];
942
943         if (j < *m) {
944
945 /*           Compute WORK(2:M) = T(J, J) L((J+1):M, J) */
946 /*            where A(J, J) = T(J, J) and A((J+1):M, J-1) = L((J+1):M, J) */
947
948             if (k > 1) {
949                 alpha = -a[j + k * a_dim1];
950                 i__1 = *m - j;
951                 daxpy_(&i__1, &alpha, &a[j + 1 + (k - 1) * a_dim1], &c__1, &
952                         work[2], &c__1);
953             }
954
955 /*           Find f2cmax(|WORK(2:M)|) */
956
957             i__1 = *m - j;
958             i2 = idamax_(&i__1, &work[2], &c__1) + 1;
959             piv = work[i2];
960
961 /*           Apply symmetric pivot */
962
963             if (i2 != 2 && piv != 0.) {
964
965 /*              Swap WORK(I1) and WORK(I2) */
966
967                 i1 = 2;
968                 work[i2] = work[i1];
969                 work[i1] = piv;
970
971 /*              Swap A(I1+1:M, I1) with A(I2, I1+1:M) */
972
973                 i1 = i1 + j - 1;
974                 i2 = i2 + j - 1;
975                 i__1 = i2 - i1 - 1;
976                 dswap_(&i__1, &a[i1 + 1 + (*j1 + i1 - 1) * a_dim1], &c__1, &a[
977                         i2 + (*j1 + i1) * a_dim1], lda);
978
979 /*              Swap A(I2+1:M, I1) with A(I2+1:M, I2) */
980
981                 if (i2 < *m) {
982                     i__1 = *m - i2;
983                     dswap_(&i__1, &a[i2 + 1 + (*j1 + i1 - 1) * a_dim1], &c__1,
984                              &a[i2 + 1 + (*j1 + i2 - 1) * a_dim1], &c__1);
985                 }
986
987 /*              Swap A(I1, I1) with A(I2, I2) */
988
989                 piv = a[i1 + (*j1 + i1 - 1) * a_dim1];
990                 a[i1 + (*j1 + i1 - 1) * a_dim1] = a[i2 + (*j1 + i2 - 1) * 
991                         a_dim1];
992                 a[i2 + (*j1 + i2 - 1) * a_dim1] = piv;
993
994 /*              Swap H(I1, I1:J1) with H(I2, I2:J1) */
995
996                 i__1 = i1 - 1;
997                 dswap_(&i__1, &h__[i1 + h_dim1], ldh, &h__[i2 + h_dim1], ldh);
998                 ipiv[i1] = i2;
999
1000                 if (i1 > k1 - 1) {
1001
1002 /*                 Swap L(1:I1-1, I1) with L(1:I1-1, I2), */
1003 /*                  skipping the first column */
1004
1005                     i__1 = i1 - k1 + 1;
1006                     dswap_(&i__1, &a[i1 + a_dim1], lda, &a[i2 + a_dim1], lda);
1007                 }
1008             } else {
1009                 ipiv[j + 1] = j + 1;
1010             }
1011
1012 /*           Set A(J+1, J) = T(J+1, J) */
1013
1014             a[j + 1 + k * a_dim1] = work[2];
1015
1016             if (j < *nb) {
1017
1018 /*              Copy A(J+1:M, J+1) into H(J+1:M, J), */
1019
1020                 i__1 = *m - j;
1021                 dcopy_(&i__1, &a[j + 1 + (k + 1) * a_dim1], &c__1, &h__[j + 1 
1022                         + (j + 1) * h_dim1], &c__1);
1023             }
1024
1025 /*           Compute L(J+2, J+1) = WORK( 3:M ) / T(J, J+1), */
1026 /*            where A(J, J+1) = T(J, J+1) and A(J+2:M, J) = L(J+2:M, J+1) */
1027
1028             if (j < *m - 1) {
1029                 if (a[j + 1 + k * a_dim1] != 0.) {
1030                     alpha = 1. / a[j + 1 + k * a_dim1];
1031                     i__1 = *m - j - 1;
1032                     dcopy_(&i__1, &work[3], &c__1, &a[j + 2 + k * a_dim1], &
1033                             c__1);
1034                     i__1 = *m - j - 1;
1035                     dscal_(&i__1, &alpha, &a[j + 2 + k * a_dim1], &c__1);
1036                 } else {
1037                     i__1 = *m - j - 1;
1038                     dlaset_("Full", &i__1, &c__1, &c_b22, &c_b22, &a[j + 2 + 
1039                             k * a_dim1], lda);
1040                 }
1041             }
1042         }
1043         ++j;
1044         goto L30;
1045 L40:
1046         ;
1047     }
1048     return 0;
1049
1050 /*     End of DLASYF_AA */
1051
1052 } /* dlasyf_aa__ */
1053