C_LAPACK: Fixes to make it compile with MSVC (#3605)
[platform/upstream/openblas.git] / lapack-netlib / SRC / ssyevd.c
1 #include <math.h>
2 #include <stdlib.h>
3 #include <string.h>
4 #include <stdio.h>
5 #include <complex.h>
6 #ifdef complex
7 #undef complex
8 #endif
9 #ifdef I
10 #undef I
11 #endif
12
13 #if defined(_WIN64)
14 typedef long long BLASLONG;
15 typedef unsigned long long BLASULONG;
16 #else
17 typedef long BLASLONG;
18 typedef unsigned long BLASULONG;
19 #endif
20
21 #ifdef LAPACK_ILP64
22 typedef BLASLONG blasint;
23 #if defined(_WIN64)
24 #define blasabs(x) llabs(x)
25 #else
26 #define blasabs(x) labs(x)
27 #endif
28 #else
29 typedef int blasint;
30 #define blasabs(x) abs(x)
31 #endif
32
33 typedef blasint integer;
34
35 typedef unsigned int uinteger;
36 typedef char *address;
37 typedef short int shortint;
38 typedef float real;
39 typedef double doublereal;
40 typedef struct { real r, i; } complex;
41 typedef struct { doublereal r, i; } doublecomplex;
42 #ifdef _MSC_VER
43 static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
44 static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
45 static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
46 static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
47 #else
48 static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
49 static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
50 static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
51 static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
52 #endif
53 #define pCf(z) (*_pCf(z))
54 #define pCd(z) (*_pCd(z))
55 typedef int logical;
56 typedef short int shortlogical;
57 typedef char logical1;
58 typedef char integer1;
59
60 #define TRUE_ (1)
61 #define FALSE_ (0)
62
63 /* Extern is for use with -E */
64 #ifndef Extern
65 #define Extern extern
66 #endif
67
68 /* I/O stuff */
69
70 typedef int flag;
71 typedef int ftnlen;
72 typedef int ftnint;
73
74 /*external read, write*/
75 typedef struct
76 {       flag cierr;
77         ftnint ciunit;
78         flag ciend;
79         char *cifmt;
80         ftnint cirec;
81 } cilist;
82
83 /*internal read, write*/
84 typedef struct
85 {       flag icierr;
86         char *iciunit;
87         flag iciend;
88         char *icifmt;
89         ftnint icirlen;
90         ftnint icirnum;
91 } icilist;
92
93 /*open*/
94 typedef struct
95 {       flag oerr;
96         ftnint ounit;
97         char *ofnm;
98         ftnlen ofnmlen;
99         char *osta;
100         char *oacc;
101         char *ofm;
102         ftnint orl;
103         char *oblnk;
104 } olist;
105
106 /*close*/
107 typedef struct
108 {       flag cerr;
109         ftnint cunit;
110         char *csta;
111 } cllist;
112
113 /*rewind, backspace, endfile*/
114 typedef struct
115 {       flag aerr;
116         ftnint aunit;
117 } alist;
118
119 /* inquire */
120 typedef struct
121 {       flag inerr;
122         ftnint inunit;
123         char *infile;
124         ftnlen infilen;
125         ftnint  *inex;  /*parameters in standard's order*/
126         ftnint  *inopen;
127         ftnint  *innum;
128         ftnint  *innamed;
129         char    *inname;
130         ftnlen  innamlen;
131         char    *inacc;
132         ftnlen  inacclen;
133         char    *inseq;
134         ftnlen  inseqlen;
135         char    *indir;
136         ftnlen  indirlen;
137         char    *infmt;
138         ftnlen  infmtlen;
139         char    *inform;
140         ftnint  informlen;
141         char    *inunf;
142         ftnlen  inunflen;
143         ftnint  *inrecl;
144         ftnint  *innrec;
145         char    *inblank;
146         ftnlen  inblanklen;
147 } inlist;
148
149 #define VOID void
150
151 union Multitype {       /* for multiple entry points */
152         integer1 g;
153         shortint h;
154         integer i;
155         /* longint j; */
156         real r;
157         doublereal d;
158         complex c;
159         doublecomplex z;
160         };
161
162 typedef union Multitype Multitype;
163
164 struct Vardesc {        /* for Namelist */
165         char *name;
166         char *addr;
167         ftnlen *dims;
168         int  type;
169         };
170 typedef struct Vardesc Vardesc;
171
172 struct Namelist {
173         char *name;
174         Vardesc **vars;
175         int nvars;
176         };
177 typedef struct Namelist Namelist;
178
179 #define abs(x) ((x) >= 0 ? (x) : -(x))
180 #define dabs(x) (fabs(x))
181 #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
182 #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
183 #define dmin(a,b) (f2cmin(a,b))
184 #define dmax(a,b) (f2cmax(a,b))
185 #define bit_test(a,b)   ((a) >> (b) & 1)
186 #define bit_clear(a,b)  ((a) & ~((uinteger)1 << (b)))
187 #define bit_set(a,b)    ((a) |  ((uinteger)1 << (b)))
188
189 #define abort_() { sig_die("Fortran abort routine called", 1); }
190 #define c_abs(z) (cabsf(Cf(z)))
191 #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
192 #ifdef _MSC_VER
193 #define c_div(c, a, b) {Cf(c)._Val[0] = (Cf(a)._Val[0]/Cf(b)._Val[0]); Cf(c)._Val[1]=(Cf(a)._Val[1]/Cf(b)._Val[1]);}
194 #define z_div(c, a, b) {Cd(c)._Val[0] = (Cd(a)._Val[0]/Cd(b)._Val[0]); Cd(c)._Val[1]=(Cd(a)._Val[1]/df(b)._Val[1]);}
195 #else
196 #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
197 #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
198 #endif
199 #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
200 #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
201 #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
202 //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
203 #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
204 #define d_abs(x) (fabs(*(x)))
205 #define d_acos(x) (acos(*(x)))
206 #define d_asin(x) (asin(*(x)))
207 #define d_atan(x) (atan(*(x)))
208 #define d_atn2(x, y) (atan2(*(x),*(y)))
209 #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
210 #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
211 #define d_cos(x) (cos(*(x)))
212 #define d_cosh(x) (cosh(*(x)))
213 #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
214 #define d_exp(x) (exp(*(x)))
215 #define d_imag(z) (cimag(Cd(z)))
216 #define r_imag(z) (cimagf(Cf(z)))
217 #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
218 #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
219 #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
220 #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
221 #define d_log(x) (log(*(x)))
222 #define d_mod(x, y) (fmod(*(x), *(y)))
223 #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
224 #define d_nint(x) u_nint(*(x))
225 #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
226 #define d_sign(a,b) u_sign(*(a),*(b))
227 #define r_sign(a,b) u_sign(*(a),*(b))
228 #define d_sin(x) (sin(*(x)))
229 #define d_sinh(x) (sinh(*(x)))
230 #define d_sqrt(x) (sqrt(*(x)))
231 #define d_tan(x) (tan(*(x)))
232 #define d_tanh(x) (tanh(*(x)))
233 #define i_abs(x) abs(*(x))
234 #define i_dnnt(x) ((integer)u_nint(*(x)))
235 #define i_len(s, n) (n)
236 #define i_nint(x) ((integer)u_nint(*(x)))
237 #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
238 #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
239 #define pow_si(B,E) spow_ui(*(B),*(E))
240 #define pow_ri(B,E) spow_ui(*(B),*(E))
241 #define pow_di(B,E) dpow_ui(*(B),*(E))
242 #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
243 #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
244 #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
245 #define s_cat(lpp, rpp, rnp, np, llp) {         ftnlen i, nc, ll; char *f__rp, *lp;     ll = (llp); lp = (lpp);         for(i=0; i < (int)*(np); ++i) {                 nc = ll;                if((rnp)[i] < nc) nc = (rnp)[i];                ll -= nc;               f__rp = (rpp)[i];               while(--nc >= 0) *lp++ = *(f__rp)++;         }  while(--ll >= 0) *lp++ = ' '; }
246 #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
247 #define s_copy(A,B,C,D) { int __i,__m; for (__i=0, __m=f2cmin((C),(D)); __i<__m && (B)[__i] != 0; ++__i) (A)[__i] = (B)[__i]; }
248 #define sig_die(s, kill) { exit(1); }
249 #define s_stop(s, n) {exit(0);}
250 static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
251 #define z_abs(z) (cabs(Cd(z)))
252 #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
253 #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
254 #define myexit_() break;
255 #define mycycle() continue;
256 #define myceiling(w) {ceil(w)}
257 #define myhuge(w) {HUGE_VAL}
258 //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
259 #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
260
261 /* procedure parameter types for -A and -C++ */
262
263 #define F2C_proc_par_types 1
264 #ifdef __cplusplus
265 typedef logical (*L_fp)(...);
266 #else
267 typedef logical (*L_fp)();
268 #endif
269
270 static float spow_ui(float x, integer n) {
271         float pow=1.0; unsigned long int u;
272         if(n != 0) {
273                 if(n < 0) n = -n, x = 1/x;
274                 for(u = n; ; ) {
275                         if(u & 01) pow *= x;
276                         if(u >>= 1) x *= x;
277                         else break;
278                 }
279         }
280         return pow;
281 }
282 static double dpow_ui(double x, integer n) {
283         double pow=1.0; unsigned long int u;
284         if(n != 0) {
285                 if(n < 0) n = -n, x = 1/x;
286                 for(u = n; ; ) {
287                         if(u & 01) pow *= x;
288                         if(u >>= 1) x *= x;
289                         else break;
290                 }
291         }
292         return pow;
293 }
294 #ifdef _MSC_VER
295 static _Fcomplex cpow_ui(complex x, integer n) {
296         complex pow={1.0,0.0}; unsigned long int u;
297                 if(n != 0) {
298                 if(n < 0) n = -n, x.r = 1/x.r, x.i=1/x.i;
299                 for(u = n; ; ) {
300                         if(u & 01) pow.r *= x.r, pow.i *= x.i;
301                         if(u >>= 1) x.r *= x.r, x.i *= x.i;
302                         else break;
303                 }
304         }
305         _Fcomplex p={pow.r, pow.i};
306         return p;
307 }
308 #else
309 static _Complex float cpow_ui(_Complex float x, integer n) {
310         _Complex float pow=1.0; unsigned long int u;
311         if(n != 0) {
312                 if(n < 0) n = -n, x = 1/x;
313                 for(u = n; ; ) {
314                         if(u & 01) pow *= x;
315                         if(u >>= 1) x *= x;
316                         else break;
317                 }
318         }
319         return pow;
320 }
321 #endif
322 #ifdef _MSC_VER
323 static _Dcomplex zpow_ui(_Dcomplex x, integer n) {
324         _Dcomplex pow={1.0,0.0}; unsigned long int u;
325         if(n != 0) {
326                 if(n < 0) n = -n, x._Val[0] = 1/x._Val[0], x._Val[1] =1/x._Val[1];
327                 for(u = n; ; ) {
328                         if(u & 01) pow._Val[0] *= x._Val[0], pow._Val[1] *= x._Val[1];
329                         if(u >>= 1) x._Val[0] *= x._Val[0], x._Val[1] *= x._Val[1];
330                         else break;
331                 }
332         }
333         _Dcomplex p = {pow._Val[0], pow._Val[1]};
334         return p;
335 }
336 #else
337 static _Complex double zpow_ui(_Complex double x, integer n) {
338         _Complex double pow=1.0; unsigned long int u;
339         if(n != 0) {
340                 if(n < 0) n = -n, x = 1/x;
341                 for(u = n; ; ) {
342                         if(u & 01) pow *= x;
343                         if(u >>= 1) x *= x;
344                         else break;
345                 }
346         }
347         return pow;
348 }
349 #endif
350 static integer pow_ii(integer x, integer n) {
351         integer pow; unsigned long int u;
352         if (n <= 0) {
353                 if (n == 0 || x == 1) pow = 1;
354                 else if (x != -1) pow = x == 0 ? 1/x : 0;
355                 else n = -n;
356         }
357         if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
358                 u = n;
359                 for(pow = 1; ; ) {
360                         if(u & 01) pow *= x;
361                         if(u >>= 1) x *= x;
362                         else break;
363                 }
364         }
365         return pow;
366 }
367 static integer dmaxloc_(double *w, integer s, integer e, integer *n)
368 {
369         double m; integer i, mi;
370         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
371                 if (w[i-1]>m) mi=i ,m=w[i-1];
372         return mi-s+1;
373 }
374 static integer smaxloc_(float *w, integer s, integer e, integer *n)
375 {
376         float m; integer i, mi;
377         for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
378                 if (w[i-1]>m) mi=i ,m=w[i-1];
379         return mi-s+1;
380 }
381 static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
382         integer n = *n_, incx = *incx_, incy = *incy_, i;
383 #ifdef _MSC_VER
384         _Fcomplex zdotc = {0.0, 0.0};
385         if (incx == 1 && incy == 1) {
386                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
387                         zdotc._Val[0] += conjf(Cf(&x[i]))._Val[0] * Cf(&y[i])._Val[0];
388                         zdotc._Val[1] += conjf(Cf(&x[i]))._Val[1] * Cf(&y[i])._Val[1];
389                 }
390         } else {
391                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
392                         zdotc._Val[0] += conjf(Cf(&x[i*incx]))._Val[0] * Cf(&y[i*incy])._Val[0];
393                         zdotc._Val[1] += conjf(Cf(&x[i*incx]))._Val[1] * Cf(&y[i*incy])._Val[1];
394                 }
395         }
396         pCf(z) = zdotc;
397 }
398 #else
399         _Complex float zdotc = 0.0;
400         if (incx == 1 && incy == 1) {
401                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
402                         zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
403                 }
404         } else {
405                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
406                         zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
407                 }
408         }
409         pCf(z) = zdotc;
410 }
411 #endif
412 static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
413         integer n = *n_, incx = *incx_, incy = *incy_, i;
414 #ifdef _MSC_VER
415         _Dcomplex zdotc = {0.0, 0.0};
416         if (incx == 1 && incy == 1) {
417                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
418                         zdotc._Val[0] += conj(Cd(&x[i]))._Val[0] * Cd(&y[i])._Val[0];
419                         zdotc._Val[1] += conj(Cd(&x[i]))._Val[1] * Cd(&y[i])._Val[1];
420                 }
421         } else {
422                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
423                         zdotc._Val[0] += conj(Cd(&x[i*incx]))._Val[0] * Cd(&y[i*incy])._Val[0];
424                         zdotc._Val[1] += conj(Cd(&x[i*incx]))._Val[1] * Cd(&y[i*incy])._Val[1];
425                 }
426         }
427         pCd(z) = zdotc;
428 }
429 #else
430         _Complex double zdotc = 0.0;
431         if (incx == 1 && incy == 1) {
432                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
433                         zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
434                 }
435         } else {
436                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
437                         zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
438                 }
439         }
440         pCd(z) = zdotc;
441 }
442 #endif  
443 static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
444         integer n = *n_, incx = *incx_, incy = *incy_, i;
445 #ifdef _MSC_VER
446         _Fcomplex zdotc = {0.0, 0.0};
447         if (incx == 1 && incy == 1) {
448                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
449                         zdotc._Val[0] += Cf(&x[i])._Val[0] * Cf(&y[i])._Val[0];
450                         zdotc._Val[1] += Cf(&x[i])._Val[1] * Cf(&y[i])._Val[1];
451                 }
452         } else {
453                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
454                         zdotc._Val[0] += Cf(&x[i*incx])._Val[0] * Cf(&y[i*incy])._Val[0];
455                         zdotc._Val[1] += Cf(&x[i*incx])._Val[1] * Cf(&y[i*incy])._Val[1];
456                 }
457         }
458         pCf(z) = zdotc;
459 }
460 #else
461         _Complex float zdotc = 0.0;
462         if (incx == 1 && incy == 1) {
463                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
464                         zdotc += Cf(&x[i]) * Cf(&y[i]);
465                 }
466         } else {
467                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
468                         zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
469                 }
470         }
471         pCf(z) = zdotc;
472 }
473 #endif
474 static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
475         integer n = *n_, incx = *incx_, incy = *incy_, i;
476 #ifdef _MSC_VER
477         _Dcomplex zdotc = {0.0, 0.0};
478         if (incx == 1 && incy == 1) {
479                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
480                         zdotc._Val[0] += Cd(&x[i])._Val[0] * Cd(&y[i])._Val[0];
481                         zdotc._Val[1] += Cd(&x[i])._Val[1] * Cd(&y[i])._Val[1];
482                 }
483         } else {
484                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
485                         zdotc._Val[0] += Cd(&x[i*incx])._Val[0] * Cd(&y[i*incy])._Val[0];
486                         zdotc._Val[1] += Cd(&x[i*incx])._Val[1] * Cd(&y[i*incy])._Val[1];
487                 }
488         }
489         pCd(z) = zdotc;
490 }
491 #else
492         _Complex double zdotc = 0.0;
493         if (incx == 1 && incy == 1) {
494                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
495                         zdotc += Cd(&x[i]) * Cd(&y[i]);
496                 }
497         } else {
498                 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
499                         zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
500                 }
501         }
502         pCd(z) = zdotc;
503 }
504 #endif
505 /*  -- translated by f2c (version 20000121).
506    You must link the resulting object file with the libraries:
507         -lf2c -lm   (in that order)
508 */
509
510
511
512
513 /* Table of constant values */
514
515 static integer c__1 = 1;
516 static integer c_n1 = -1;
517 static integer c__0 = 0;
518 static real c_b17 = 1.f;
519
520 /* > \brief <b> SSYEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY mat
521 rices</b> */
522
523 /*  =========== DOCUMENTATION =========== */
524
525 /* Online html documentation available at */
526 /*            http://www.netlib.org/lapack/explore-html/ */
527
528 /* > \htmlonly */
529 /* > Download SSYEVD + dependencies */
530 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ssyevd.
531 f"> */
532 /* > [TGZ]</a> */
533 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ssyevd.
534 f"> */
535 /* > [ZIP]</a> */
536 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ssyevd.
537 f"> */
538 /* > [TXT]</a> */
539 /* > \endhtmlonly */
540
541 /*  Definition: */
542 /*  =========== */
543
544 /*       SUBROUTINE SSYEVD( JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, IWORK, */
545 /*                          LIWORK, INFO ) */
546
547 /*       CHARACTER          JOBZ, UPLO */
548 /*       INTEGER            INFO, LDA, LIWORK, LWORK, N */
549 /*       INTEGER            IWORK( * ) */
550 /*       REAL               A( LDA, * ), W( * ), WORK( * ) */
551
552
553 /* > \par Purpose: */
554 /*  ============= */
555 /* > */
556 /* > \verbatim */
557 /* > */
558 /* > SSYEVD computes all eigenvalues and, optionally, eigenvectors of a */
559 /* > real symmetric matrix A. If eigenvectors are desired, it uses a */
560 /* > divide and conquer algorithm. */
561 /* > */
562 /* > The divide and conquer algorithm makes very mild assumptions about */
563 /* > floating point arithmetic. It will work on machines with a guard */
564 /* > digit in add/subtract, or on those binary machines without guard */
565 /* > digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or */
566 /* > Cray-2. It could conceivably fail on hexadecimal or decimal machines */
567 /* > without guard digits, but we know of none. */
568 /* > */
569 /* > Because of large use of BLAS of level 3, SSYEVD needs N**2 more */
570 /* > workspace than SSYEVX. */
571 /* > \endverbatim */
572
573 /*  Arguments: */
574 /*  ========== */
575
576 /* > \param[in] JOBZ */
577 /* > \verbatim */
578 /* >          JOBZ is CHARACTER*1 */
579 /* >          = 'N':  Compute eigenvalues only; */
580 /* >          = 'V':  Compute eigenvalues and eigenvectors. */
581 /* > \endverbatim */
582 /* > */
583 /* > \param[in] UPLO */
584 /* > \verbatim */
585 /* >          UPLO is CHARACTER*1 */
586 /* >          = 'U':  Upper triangle of A is stored; */
587 /* >          = 'L':  Lower triangle of A is stored. */
588 /* > \endverbatim */
589 /* > */
590 /* > \param[in] N */
591 /* > \verbatim */
592 /* >          N is INTEGER */
593 /* >          The order of the matrix A.  N >= 0. */
594 /* > \endverbatim */
595 /* > */
596 /* > \param[in,out] A */
597 /* > \verbatim */
598 /* >          A is REAL array, dimension (LDA, N) */
599 /* >          On entry, the symmetric matrix A.  If UPLO = 'U', the */
600 /* >          leading N-by-N upper triangular part of A contains the */
601 /* >          upper triangular part of the matrix A.  If UPLO = 'L', */
602 /* >          the leading N-by-N lower triangular part of A contains */
603 /* >          the lower triangular part of the matrix A. */
604 /* >          On exit, if JOBZ = 'V', then if INFO = 0, A contains the */
605 /* >          orthonormal eigenvectors of the matrix A. */
606 /* >          If JOBZ = 'N', then on exit the lower triangle (if UPLO='L') */
607 /* >          or the upper triangle (if UPLO='U') of A, including the */
608 /* >          diagonal, is destroyed. */
609 /* > \endverbatim */
610 /* > */
611 /* > \param[in] LDA */
612 /* > \verbatim */
613 /* >          LDA is INTEGER */
614 /* >          The leading dimension of the array A.  LDA >= f2cmax(1,N). */
615 /* > \endverbatim */
616 /* > */
617 /* > \param[out] W */
618 /* > \verbatim */
619 /* >          W is REAL array, dimension (N) */
620 /* >          If INFO = 0, the eigenvalues in ascending order. */
621 /* > \endverbatim */
622 /* > */
623 /* > \param[out] WORK */
624 /* > \verbatim */
625 /* >          WORK is REAL array, */
626 /* >                                         dimension (LWORK) */
627 /* >          On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
628 /* > \endverbatim */
629 /* > */
630 /* > \param[in] LWORK */
631 /* > \verbatim */
632 /* >          LWORK is INTEGER */
633 /* >          The dimension of the array WORK. */
634 /* >          If N <= 1,               LWORK must be at least 1. */
635 /* >          If JOBZ = 'N' and N > 1, LWORK must be at least 2*N+1. */
636 /* >          If JOBZ = 'V' and N > 1, LWORK must be at least */
637 /* >                                                1 + 6*N + 2*N**2. */
638 /* > */
639 /* >          If LWORK = -1, then a workspace query is assumed; the routine */
640 /* >          only calculates the optimal sizes of the WORK and IWORK */
641 /* >          arrays, returns these values as the first entries of the WORK */
642 /* >          and IWORK arrays, and no error message related to LWORK or */
643 /* >          LIWORK is issued by XERBLA. */
644 /* > \endverbatim */
645 /* > */
646 /* > \param[out] IWORK */
647 /* > \verbatim */
648 /* >          IWORK is INTEGER array, dimension (MAX(1,LIWORK)) */
649 /* >          On exit, if INFO = 0, IWORK(1) returns the optimal LIWORK. */
650 /* > \endverbatim */
651 /* > */
652 /* > \param[in] LIWORK */
653 /* > \verbatim */
654 /* >          LIWORK is INTEGER */
655 /* >          The dimension of the array IWORK. */
656 /* >          If N <= 1,                LIWORK must be at least 1. */
657 /* >          If JOBZ  = 'N' and N > 1, LIWORK must be at least 1. */
658 /* >          If JOBZ  = 'V' and N > 1, LIWORK must be at least 3 + 5*N. */
659 /* > */
660 /* >          If LIWORK = -1, then a workspace query is assumed; the */
661 /* >          routine only calculates the optimal sizes of the WORK and */
662 /* >          IWORK arrays, returns these values as the first entries of */
663 /* >          the WORK and IWORK arrays, and no error message related to */
664 /* >          LWORK or LIWORK is issued by XERBLA. */
665 /* > \endverbatim */
666 /* > */
667 /* > \param[out] INFO */
668 /* > \verbatim */
669 /* >          INFO is INTEGER */
670 /* >          = 0:  successful exit */
671 /* >          < 0:  if INFO = -i, the i-th argument had an illegal value */
672 /* >          > 0:  if INFO = i and JOBZ = 'N', then the algorithm failed */
673 /* >                to converge; i off-diagonal elements of an intermediate */
674 /* >                tridiagonal form did not converge to zero; */
675 /* >                if INFO = i and JOBZ = 'V', then the algorithm failed */
676 /* >                to compute an eigenvalue while working on the submatrix */
677 /* >                lying in rows and columns INFO/(N+1) through */
678 /* >                mod(INFO,N+1). */
679 /* > \endverbatim */
680
681 /*  Authors: */
682 /*  ======== */
683
684 /* > \author Univ. of Tennessee */
685 /* > \author Univ. of California Berkeley */
686 /* > \author Univ. of Colorado Denver */
687 /* > \author NAG Ltd. */
688
689 /* > \date December 2016 */
690
691 /* > \ingroup realSYeigen */
692
693 /* > \par Contributors: */
694 /*  ================== */
695 /* > */
696 /* > Jeff Rutter, Computer Science Division, University of California */
697 /* > at Berkeley, USA \n */
698 /* >  Modified by Francoise Tisseur, University of Tennessee \n */
699 /* >  Modified description of INFO. Sven, 16 Feb 05. \n */
700 /* > */
701 /*  ===================================================================== */
702 /* Subroutine */ int ssyevd_(char *jobz, char *uplo, integer *n, real *a, 
703         integer *lda, real *w, real *work, integer *lwork, integer *iwork, 
704         integer *liwork, integer *info)
705 {
706     /* System generated locals */
707     integer a_dim1, a_offset, i__1, i__2;
708     real r__1;
709
710     /* Local variables */
711     integer inde;
712     real anrm, rmin, rmax;
713     integer lopt;
714     real sigma;
715     extern logical lsame_(char *, char *);
716     integer iinfo;
717     extern /* Subroutine */ int sscal_(integer *, real *, real *, integer *);
718     integer lwmin, liopt;
719     logical lower, wantz;
720     integer indwk2, llwrk2, iscale;
721     extern real slamch_(char *);
722     real safmin;
723     extern integer ilaenv_(integer *, char *, char *, integer *, integer *, 
724             integer *, integer *, ftnlen, ftnlen);
725     extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
726     real bignum;
727     extern /* Subroutine */ int slascl_(char *, integer *, integer *, real *, 
728             real *, integer *, integer *, real *, integer *, integer *);
729     integer indtau;
730     extern /* Subroutine */ int sstedc_(char *, integer *, real *, real *, 
731             real *, integer *, real *, integer *, integer *, integer *, 
732             integer *), slacpy_(char *, integer *, integer *, real *, 
733             integer *, real *, integer *);
734     integer indwrk, liwmin;
735     extern /* Subroutine */ int ssterf_(integer *, real *, real *, integer *);
736     extern real slansy_(char *, char *, integer *, real *, integer *, real *);
737     integer llwork;
738     real smlnum;
739     logical lquery;
740     extern /* Subroutine */ int sormtr_(char *, char *, char *, integer *, 
741             integer *, real *, integer *, real *, real *, integer *, real *, 
742             integer *, integer *), ssytrd_(char *, 
743             integer *, real *, integer *, real *, real *, real *, real *, 
744             integer *, integer *);
745     real eps;
746
747
748 /*  -- LAPACK driver routine (version 3.7.0) -- */
749 /*  -- LAPACK is a software package provided by Univ. of Tennessee,    -- */
750 /*  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
751 /*     December 2016 */
752
753
754 /*  ===================================================================== */
755
756
757
758 /*     Test the input parameters. */
759
760     /* Parameter adjustments */
761     a_dim1 = *lda;
762     a_offset = 1 + a_dim1 * 1;
763     a -= a_offset;
764     --w;
765     --work;
766     --iwork;
767
768     /* Function Body */
769     wantz = lsame_(jobz, "V");
770     lower = lsame_(uplo, "L");
771     lquery = *lwork == -1 || *liwork == -1;
772
773     *info = 0;
774     if (! (wantz || lsame_(jobz, "N"))) {
775         *info = -1;
776     } else if (! (lower || lsame_(uplo, "U"))) {
777         *info = -2;
778     } else if (*n < 0) {
779         *info = -3;
780     } else if (*lda < f2cmax(1,*n)) {
781         *info = -5;
782     }
783
784     if (*info == 0) {
785         if (*n <= 1) {
786             liwmin = 1;
787             lwmin = 1;
788             lopt = lwmin;
789             liopt = liwmin;
790         } else {
791             if (wantz) {
792                 liwmin = *n * 5 + 3;
793 /* Computing 2nd power */
794                 i__1 = *n;
795                 lwmin = *n * 6 + 1 + (i__1 * i__1 << 1);
796             } else {
797                 liwmin = 1;
798                 lwmin = (*n << 1) + 1;
799             }
800 /* Computing MAX */
801             i__1 = lwmin, i__2 = (*n << 1) + ilaenv_(&c__1, "SSYTRD", uplo, n,
802                      &c_n1, &c_n1, &c_n1, (ftnlen)6, (ftnlen)1);
803             lopt = f2cmax(i__1,i__2);
804             liopt = liwmin;
805         }
806         work[1] = (real) lopt;
807         iwork[1] = liopt;
808
809         if (*lwork < lwmin && ! lquery) {
810             *info = -8;
811         } else if (*liwork < liwmin && ! lquery) {
812             *info = -10;
813         }
814     }
815
816     if (*info != 0) {
817         i__1 = -(*info);
818         xerbla_("SSYEVD", &i__1, (ftnlen)6);
819         return 0;
820     } else if (lquery) {
821         return 0;
822     }
823
824 /*     Quick return if possible */
825
826     if (*n == 0) {
827         return 0;
828     }
829
830     if (*n == 1) {
831         w[1] = a[a_dim1 + 1];
832         if (wantz) {
833             a[a_dim1 + 1] = 1.f;
834         }
835         return 0;
836     }
837
838 /*     Get machine constants. */
839
840     safmin = slamch_("Safe minimum");
841     eps = slamch_("Precision");
842     smlnum = safmin / eps;
843     bignum = 1.f / smlnum;
844     rmin = sqrt(smlnum);
845     rmax = sqrt(bignum);
846
847 /*     Scale matrix to allowable range, if necessary. */
848
849     anrm = slansy_("M", uplo, n, &a[a_offset], lda, &work[1]);
850     iscale = 0;
851     if (anrm > 0.f && anrm < rmin) {
852         iscale = 1;
853         sigma = rmin / anrm;
854     } else if (anrm > rmax) {
855         iscale = 1;
856         sigma = rmax / anrm;
857     }
858     if (iscale == 1) {
859         slascl_(uplo, &c__0, &c__0, &c_b17, &sigma, n, n, &a[a_offset], lda, 
860                 info);
861     }
862
863 /*     Call SSYTRD to reduce symmetric matrix to tridiagonal form. */
864
865     inde = 1;
866     indtau = inde + *n;
867     indwrk = indtau + *n;
868     llwork = *lwork - indwrk + 1;
869     indwk2 = indwrk + *n * *n;
870     llwrk2 = *lwork - indwk2 + 1;
871
872     ssytrd_(uplo, n, &a[a_offset], lda, &w[1], &work[inde], &work[indtau], &
873             work[indwrk], &llwork, &iinfo);
874
875 /*     For eigenvalues only, call SSTERF.  For eigenvectors, first call */
876 /*     SSTEDC to generate the eigenvector matrix, WORK(INDWRK), of the */
877 /*     tridiagonal matrix, then call SORMTR to multiply it by the */
878 /*     Householder transformations stored in A. */
879
880     if (! wantz) {
881         ssterf_(n, &w[1], &work[inde], info);
882     } else {
883         sstedc_("I", n, &w[1], &work[inde], &work[indwrk], n, &work[indwk2], &
884                 llwrk2, &iwork[1], liwork, info);
885         sormtr_("L", uplo, "N", n, n, &a[a_offset], lda, &work[indtau], &work[
886                 indwrk], n, &work[indwk2], &llwrk2, &iinfo);
887         slacpy_("A", n, n, &work[indwrk], n, &a[a_offset], lda);
888     }
889
890 /*     If matrix was scaled, then rescale eigenvalues appropriately. */
891
892     if (iscale == 1) {
893         r__1 = 1.f / sigma;
894         sscal_(n, &r__1, &w[1], &c__1);
895     }
896
897     work[1] = (real) lopt;
898     iwork[1] = liopt;
899
900     return 0;
901
902 /*     End of SSYEVD */
903
904 } /* ssyevd_ */
905