You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

dsyev.c 20 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715
  1. /* f2c.h -- Standard Fortran to C header file */
  2. /** barf [ba:rf] 2. "He suggested using FORTRAN, and everybody barfed."
  3. - From The Shogakukan DICTIONARY OF NEW ENGLISH (Second edition) */
  4. #ifndef F2C_INCLUDE
  5. #define F2C_INCLUDE
  6. #include <math.h>
  7. #include <stdlib.h>
  8. #include <string.h>
  9. #include <stdio.h>
  10. #include <complex.h>
  11. #ifdef complex
  12. #undef complex
  13. #endif
  14. #ifdef I
  15. #undef I
  16. #endif
  17. typedef int integer;
  18. typedef unsigned int uinteger;
  19. typedef char *address;
  20. typedef short int shortint;
  21. typedef float real;
  22. typedef double doublereal;
  23. typedef struct { real r, i; } complex;
  24. typedef struct { doublereal r, i; } doublecomplex;
  25. static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
  26. static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
  27. static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
  28. static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
  29. #define pCf(z) (*_pCf(z))
  30. #define pCd(z) (*_pCd(z))
  31. typedef int logical;
  32. typedef short int shortlogical;
  33. typedef char logical1;
  34. typedef char integer1;
  35. #define TRUE_ (1)
  36. #define FALSE_ (0)
  37. /* Extern is for use with -E */
  38. #ifndef Extern
  39. #define Extern extern
  40. #endif
  41. /* I/O stuff */
  42. typedef int flag;
  43. typedef int ftnlen;
  44. typedef int ftnint;
  45. /*external read, write*/
  46. typedef struct
  47. { flag cierr;
  48. ftnint ciunit;
  49. flag ciend;
  50. char *cifmt;
  51. ftnint cirec;
  52. } cilist;
  53. /*internal read, write*/
  54. typedef struct
  55. { flag icierr;
  56. char *iciunit;
  57. flag iciend;
  58. char *icifmt;
  59. ftnint icirlen;
  60. ftnint icirnum;
  61. } icilist;
  62. /*open*/
  63. typedef struct
  64. { flag oerr;
  65. ftnint ounit;
  66. char *ofnm;
  67. ftnlen ofnmlen;
  68. char *osta;
  69. char *oacc;
  70. char *ofm;
  71. ftnint orl;
  72. char *oblnk;
  73. } olist;
  74. /*close*/
  75. typedef struct
  76. { flag cerr;
  77. ftnint cunit;
  78. char *csta;
  79. } cllist;
  80. /*rewind, backspace, endfile*/
  81. typedef struct
  82. { flag aerr;
  83. ftnint aunit;
  84. } alist;
  85. /* inquire */
  86. typedef struct
  87. { flag inerr;
  88. ftnint inunit;
  89. char *infile;
  90. ftnlen infilen;
  91. ftnint *inex; /*parameters in standard's order*/
  92. ftnint *inopen;
  93. ftnint *innum;
  94. ftnint *innamed;
  95. char *inname;
  96. ftnlen innamlen;
  97. char *inacc;
  98. ftnlen inacclen;
  99. char *inseq;
  100. ftnlen inseqlen;
  101. char *indir;
  102. ftnlen indirlen;
  103. char *infmt;
  104. ftnlen infmtlen;
  105. char *inform;
  106. ftnint informlen;
  107. char *inunf;
  108. ftnlen inunflen;
  109. ftnint *inrecl;
  110. ftnint *innrec;
  111. char *inblank;
  112. ftnlen inblanklen;
  113. } inlist;
  114. #define VOID void
  115. union Multitype { /* for multiple entry points */
  116. integer1 g;
  117. shortint h;
  118. integer i;
  119. /* longint j; */
  120. real r;
  121. doublereal d;
  122. complex c;
  123. doublecomplex z;
  124. };
  125. typedef union Multitype Multitype;
  126. struct Vardesc { /* for Namelist */
  127. char *name;
  128. char *addr;
  129. ftnlen *dims;
  130. int type;
  131. };
  132. typedef struct Vardesc Vardesc;
  133. struct Namelist {
  134. char *name;
  135. Vardesc **vars;
  136. int nvars;
  137. };
  138. typedef struct Namelist Namelist;
  139. #define abs(x) ((x) >= 0 ? (x) : -(x))
  140. #define dabs(x) (fabs(x))
  141. #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
  142. #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
  143. #define dmin(a,b) (f2cmin(a,b))
  144. #define dmax(a,b) (f2cmax(a,b))
  145. #define bit_test(a,b) ((a) >> (b) & 1)
  146. #define bit_clear(a,b) ((a) & ~((uinteger)1 << (b)))
  147. #define bit_set(a,b) ((a) | ((uinteger)1 << (b)))
  148. #define abort_() { sig_die("Fortran abort routine called", 1); }
  149. #define c_abs(z) (cabsf(Cf(z)))
  150. #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
  151. #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
  152. #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
  153. #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
  154. #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
  155. #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
  156. //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
  157. #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
  158. #define d_abs(x) (fabs(*(x)))
  159. #define d_acos(x) (acos(*(x)))
  160. #define d_asin(x) (asin(*(x)))
  161. #define d_atan(x) (atan(*(x)))
  162. #define d_atn2(x, y) (atan2(*(x),*(y)))
  163. #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
  164. #define r_cnjg(R, Z) { pCf(R) = conj(Cf(Z)); }
  165. #define d_cos(x) (cos(*(x)))
  166. #define d_cosh(x) (cosh(*(x)))
  167. #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
  168. #define d_exp(x) (exp(*(x)))
  169. #define d_imag(z) (cimag(Cd(z)))
  170. #define r_imag(z) (cimag(Cf(z)))
  171. #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
  172. #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
  173. #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
  174. #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
  175. #define d_log(x) (log(*(x)))
  176. #define d_mod(x, y) (fmod(*(x), *(y)))
  177. #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
  178. #define d_nint(x) u_nint(*(x))
  179. #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
  180. #define d_sign(a,b) u_sign(*(a),*(b))
  181. #define r_sign(a,b) u_sign(*(a),*(b))
  182. #define d_sin(x) (sin(*(x)))
  183. #define d_sinh(x) (sinh(*(x)))
  184. #define d_sqrt(x) (sqrt(*(x)))
  185. #define d_tan(x) (tan(*(x)))
  186. #define d_tanh(x) (tanh(*(x)))
  187. #define i_abs(x) abs(*(x))
  188. #define i_dnnt(x) ((integer)u_nint(*(x)))
  189. #define i_len(s, n) (n)
  190. #define i_nint(x) ((integer)u_nint(*(x)))
  191. #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
  192. #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
  193. #define pow_si(B,E) spow_ui(*(B),*(E))
  194. #define pow_ri(B,E) spow_ui(*(B),*(E))
  195. #define pow_di(B,E) dpow_ui(*(B),*(E))
  196. #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
  197. #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
  198. #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
  199. #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++ = ' '; }
  200. #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
  201. #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]; }
  202. #define sig_die(s, kill) { exit(1); }
  203. #define s_stop(s, n) {exit(0);}
  204. static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
  205. #define z_abs(z) (cabs(Cd(z)))
  206. #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
  207. #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
  208. #define myexit_() break;
  209. #define mycycle() continue;
  210. #define myceiling(w) {ceil(w)}
  211. #define myhuge(w) {HUGE_VAL}
  212. //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
  213. #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
  214. /* procedure parameter types for -A and -C++ */
  215. #define F2C_proc_par_types 1
  216. #ifdef __cplusplus
  217. typedef logical (*L_fp)(...);
  218. #else
  219. typedef logical (*L_fp)();
  220. #endif
  221. static float spow_ui(float x, integer n) {
  222. float pow=1.0; unsigned long int u;
  223. if(n != 0) {
  224. if(n < 0) n = -n, x = 1/x;
  225. for(u = n; ; ) {
  226. if(u & 01) pow *= x;
  227. if(u >>= 1) x *= x;
  228. else break;
  229. }
  230. }
  231. return pow;
  232. }
  233. static double dpow_ui(double x, integer n) {
  234. double pow=1.0; unsigned long int u;
  235. if(n != 0) {
  236. if(n < 0) n = -n, x = 1/x;
  237. for(u = n; ; ) {
  238. if(u & 01) pow *= x;
  239. if(u >>= 1) x *= x;
  240. else break;
  241. }
  242. }
  243. return pow;
  244. }
  245. static _Complex float cpow_ui(_Complex float x, integer n) {
  246. _Complex float pow=1.0; unsigned long int u;
  247. if(n != 0) {
  248. if(n < 0) n = -n, x = 1/x;
  249. for(u = n; ; ) {
  250. if(u & 01) pow *= x;
  251. if(u >>= 1) x *= x;
  252. else break;
  253. }
  254. }
  255. return pow;
  256. }
  257. static _Complex double zpow_ui(_Complex double x, integer n) {
  258. _Complex double pow=1.0; unsigned long int u;
  259. if(n != 0) {
  260. if(n < 0) n = -n, x = 1/x;
  261. for(u = n; ; ) {
  262. if(u & 01) pow *= x;
  263. if(u >>= 1) x *= x;
  264. else break;
  265. }
  266. }
  267. return pow;
  268. }
  269. static integer pow_ii(integer x, integer n) {
  270. integer pow; unsigned long int u;
  271. if (n <= 0) {
  272. if (n == 0 || x == 1) pow = 1;
  273. else if (x != -1) pow = x == 0 ? 1/x : 0;
  274. else n = -n;
  275. }
  276. if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
  277. u = n;
  278. for(pow = 1; ; ) {
  279. if(u & 01) pow *= x;
  280. if(u >>= 1) x *= x;
  281. else break;
  282. }
  283. }
  284. return pow;
  285. }
  286. static integer dmaxloc_(double *w, integer s, integer e, integer *n)
  287. {
  288. double m; integer i, mi;
  289. for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
  290. if (w[i-1]>m) mi=i ,m=w[i-1];
  291. return mi-s+1;
  292. }
  293. static integer smaxloc_(float *w, integer s, integer e, integer *n)
  294. {
  295. float m; integer i, mi;
  296. for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
  297. if (w[i-1]>m) mi=i ,m=w[i-1];
  298. return mi-s+1;
  299. }
  300. static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
  301. integer n = *n_, incx = *incx_, incy = *incy_, i;
  302. _Complex float zdotc = 0.0;
  303. if (incx == 1 && incy == 1) {
  304. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  305. zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
  306. }
  307. } else {
  308. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  309. zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
  310. }
  311. }
  312. pCf(z) = zdotc;
  313. }
  314. static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
  315. integer n = *n_, incx = *incx_, incy = *incy_, i;
  316. _Complex double zdotc = 0.0;
  317. if (incx == 1 && incy == 1) {
  318. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  319. zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
  320. }
  321. } else {
  322. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  323. zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
  324. }
  325. }
  326. pCd(z) = zdotc;
  327. }
  328. static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
  329. integer n = *n_, incx = *incx_, incy = *incy_, i;
  330. _Complex float zdotc = 0.0;
  331. if (incx == 1 && incy == 1) {
  332. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  333. zdotc += Cf(&x[i]) * Cf(&y[i]);
  334. }
  335. } else {
  336. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  337. zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
  338. }
  339. }
  340. pCf(z) = zdotc;
  341. }
  342. static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
  343. integer n = *n_, incx = *incx_, incy = *incy_, i;
  344. _Complex double zdotc = 0.0;
  345. if (incx == 1 && incy == 1) {
  346. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  347. zdotc += Cd(&x[i]) * Cd(&y[i]);
  348. }
  349. } else {
  350. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  351. zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
  352. }
  353. }
  354. pCd(z) = zdotc;
  355. }
  356. #endif
  357. /* -- translated by f2c (version 20000121).
  358. You must link the resulting object file with the libraries:
  359. -lf2c -lm (in that order)
  360. */
  361. /* Table of constant values */
  362. static integer c__1 = 1;
  363. static integer c_n1 = -1;
  364. static integer c__0 = 0;
  365. static doublereal c_b17 = 1.;
  366. /* > \brief <b> DSYEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for SY matr
  367. ices</b> */
  368. /* =========== DOCUMENTATION =========== */
  369. /* Online html documentation available at */
  370. /* http://www.netlib.org/lapack/explore-html/ */
  371. /* > \htmlonly */
  372. /* > Download DSYEV + dependencies */
  373. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsyev.f
  374. "> */
  375. /* > [TGZ]</a> */
  376. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsyev.f
  377. "> */
  378. /* > [ZIP]</a> */
  379. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsyev.f
  380. "> */
  381. /* > [TXT]</a> */
  382. /* > \endhtmlonly */
  383. /* Definition: */
  384. /* =========== */
  385. /* SUBROUTINE DSYEV( JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, INFO ) */
  386. /* CHARACTER JOBZ, UPLO */
  387. /* INTEGER INFO, LDA, LWORK, N */
  388. /* DOUBLE PRECISION A( LDA, * ), W( * ), WORK( * ) */
  389. /* > \par Purpose: */
  390. /* ============= */
  391. /* > */
  392. /* > \verbatim */
  393. /* > */
  394. /* > DSYEV computes all eigenvalues and, optionally, eigenvectors of a */
  395. /* > real symmetric matrix A. */
  396. /* > \endverbatim */
  397. /* Arguments: */
  398. /* ========== */
  399. /* > \param[in] JOBZ */
  400. /* > \verbatim */
  401. /* > JOBZ is CHARACTER*1 */
  402. /* > = 'N': Compute eigenvalues only; */
  403. /* > = 'V': Compute eigenvalues and eigenvectors. */
  404. /* > \endverbatim */
  405. /* > */
  406. /* > \param[in] UPLO */
  407. /* > \verbatim */
  408. /* > UPLO is CHARACTER*1 */
  409. /* > = 'U': Upper triangle of A is stored; */
  410. /* > = 'L': Lower triangle of A is stored. */
  411. /* > \endverbatim */
  412. /* > */
  413. /* > \param[in] N */
  414. /* > \verbatim */
  415. /* > N is INTEGER */
  416. /* > The order of the matrix A. N >= 0. */
  417. /* > \endverbatim */
  418. /* > */
  419. /* > \param[in,out] A */
  420. /* > \verbatim */
  421. /* > A is DOUBLE PRECISION array, dimension (LDA, N) */
  422. /* > On entry, the symmetric matrix A. If UPLO = 'U', the */
  423. /* > leading N-by-N upper triangular part of A contains the */
  424. /* > upper triangular part of the matrix A. If UPLO = 'L', */
  425. /* > the leading N-by-N lower triangular part of A contains */
  426. /* > the lower triangular part of the matrix A. */
  427. /* > On exit, if JOBZ = 'V', then if INFO = 0, A contains the */
  428. /* > orthonormal eigenvectors of the matrix A. */
  429. /* > If JOBZ = 'N', then on exit the lower triangle (if UPLO='L') */
  430. /* > or the upper triangle (if UPLO='U') of A, including the */
  431. /* > diagonal, is destroyed. */
  432. /* > \endverbatim */
  433. /* > */
  434. /* > \param[in] LDA */
  435. /* > \verbatim */
  436. /* > LDA is INTEGER */
  437. /* > The leading dimension of the array A. LDA >= f2cmax(1,N). */
  438. /* > \endverbatim */
  439. /* > */
  440. /* > \param[out] W */
  441. /* > \verbatim */
  442. /* > W is DOUBLE PRECISION array, dimension (N) */
  443. /* > If INFO = 0, the eigenvalues in ascending order. */
  444. /* > \endverbatim */
  445. /* > */
  446. /* > \param[out] WORK */
  447. /* > \verbatim */
  448. /* > WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK)) */
  449. /* > On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
  450. /* > \endverbatim */
  451. /* > */
  452. /* > \param[in] LWORK */
  453. /* > \verbatim */
  454. /* > LWORK is INTEGER */
  455. /* > The length of the array WORK. LWORK >= f2cmax(1,3*N-1). */
  456. /* > For optimal efficiency, LWORK >= (NB+2)*N, */
  457. /* > where NB is the blocksize for DSYTRD returned by ILAENV. */
  458. /* > */
  459. /* > If LWORK = -1, then a workspace query is assumed; the routine */
  460. /* > only calculates the optimal size of the WORK array, returns */
  461. /* > this value as the first entry of the WORK array, and no error */
  462. /* > message related to LWORK is issued by XERBLA. */
  463. /* > \endverbatim */
  464. /* > */
  465. /* > \param[out] INFO */
  466. /* > \verbatim */
  467. /* > INFO is INTEGER */
  468. /* > = 0: successful exit */
  469. /* > < 0: if INFO = -i, the i-th argument had an illegal value */
  470. /* > > 0: if INFO = i, the algorithm failed to converge; i */
  471. /* > off-diagonal elements of an intermediate tridiagonal */
  472. /* > form did not converge to zero. */
  473. /* > \endverbatim */
  474. /* Authors: */
  475. /* ======== */
  476. /* > \author Univ. of Tennessee */
  477. /* > \author Univ. of California Berkeley */
  478. /* > \author Univ. of Colorado Denver */
  479. /* > \author NAG Ltd. */
  480. /* > \date December 2016 */
  481. /* > \ingroup doubleSYeigen */
  482. /* ===================================================================== */
  483. /* Subroutine */ int dsyev_(char *jobz, char *uplo, integer *n, doublereal *a,
  484. integer *lda, doublereal *w, doublereal *work, integer *lwork,
  485. integer *info)
  486. {
  487. /* System generated locals */
  488. integer a_dim1, a_offset, i__1, i__2;
  489. doublereal d__1;
  490. /* Local variables */
  491. integer inde;
  492. doublereal anrm;
  493. integer imax;
  494. doublereal rmin, rmax;
  495. extern /* Subroutine */ int dscal_(integer *, doublereal *, doublereal *,
  496. integer *);
  497. doublereal sigma;
  498. extern logical lsame_(char *, char *);
  499. integer iinfo;
  500. logical lower, wantz;
  501. integer nb;
  502. extern doublereal dlamch_(char *);
  503. integer iscale;
  504. extern /* Subroutine */ int dlascl_(char *, integer *, integer *,
  505. doublereal *, doublereal *, integer *, integer *, doublereal *,
  506. integer *, integer *);
  507. doublereal safmin;
  508. extern integer ilaenv_(integer *, char *, char *, integer *, integer *,
  509. integer *, integer *, ftnlen, ftnlen);
  510. extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
  511. doublereal bignum;
  512. integer indtau;
  513. extern /* Subroutine */ int dsterf_(integer *, doublereal *, doublereal *,
  514. integer *);
  515. extern doublereal dlansy_(char *, char *, integer *, doublereal *,
  516. integer *, doublereal *);
  517. integer indwrk;
  518. extern /* Subroutine */ int dorgtr_(char *, integer *, doublereal *,
  519. integer *, doublereal *, doublereal *, integer *, integer *), dsteqr_(char *, integer *, doublereal *, doublereal *,
  520. doublereal *, integer *, doublereal *, integer *),
  521. dsytrd_(char *, integer *, doublereal *, integer *, doublereal *,
  522. doublereal *, doublereal *, doublereal *, integer *, integer *);
  523. integer llwork;
  524. doublereal smlnum;
  525. integer lwkopt;
  526. logical lquery;
  527. doublereal eps;
  528. /* -- LAPACK driver routine (version 3.7.0) -- */
  529. /* -- LAPACK is a software package provided by Univ. of Tennessee, -- */
  530. /* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
  531. /* December 2016 */
  532. /* ===================================================================== */
  533. /* Test the input parameters. */
  534. /* Parameter adjustments */
  535. a_dim1 = *lda;
  536. a_offset = 1 + a_dim1 * 1;
  537. a -= a_offset;
  538. --w;
  539. --work;
  540. /* Function Body */
  541. wantz = lsame_(jobz, "V");
  542. lower = lsame_(uplo, "L");
  543. lquery = *lwork == -1;
  544. *info = 0;
  545. if (! (wantz || lsame_(jobz, "N"))) {
  546. *info = -1;
  547. } else if (! (lower || lsame_(uplo, "U"))) {
  548. *info = -2;
  549. } else if (*n < 0) {
  550. *info = -3;
  551. } else if (*lda < f2cmax(1,*n)) {
  552. *info = -5;
  553. }
  554. if (*info == 0) {
  555. nb = ilaenv_(&c__1, "DSYTRD", uplo, n, &c_n1, &c_n1, &c_n1, (ftnlen)6,
  556. (ftnlen)1);
  557. /* Computing MAX */
  558. i__1 = 1, i__2 = (nb + 2) * *n;
  559. lwkopt = f2cmax(i__1,i__2);
  560. work[1] = (doublereal) lwkopt;
  561. /* Computing MAX */
  562. i__1 = 1, i__2 = *n * 3 - 1;
  563. if (*lwork < f2cmax(i__1,i__2) && ! lquery) {
  564. *info = -8;
  565. }
  566. }
  567. if (*info != 0) {
  568. i__1 = -(*info);
  569. xerbla_("DSYEV ", &i__1, (ftnlen)6);
  570. return 0;
  571. } else if (lquery) {
  572. return 0;
  573. }
  574. /* Quick return if possible */
  575. if (*n == 0) {
  576. return 0;
  577. }
  578. if (*n == 1) {
  579. w[1] = a[a_dim1 + 1];
  580. work[1] = 2.;
  581. if (wantz) {
  582. a[a_dim1 + 1] = 1.;
  583. }
  584. return 0;
  585. }
  586. /* Get machine constants. */
  587. safmin = dlamch_("Safe minimum");
  588. eps = dlamch_("Precision");
  589. smlnum = safmin / eps;
  590. bignum = 1. / smlnum;
  591. rmin = sqrt(smlnum);
  592. rmax = sqrt(bignum);
  593. /* Scale matrix to allowable range, if necessary. */
  594. anrm = dlansy_("M", uplo, n, &a[a_offset], lda, &work[1]);
  595. iscale = 0;
  596. if (anrm > 0. && anrm < rmin) {
  597. iscale = 1;
  598. sigma = rmin / anrm;
  599. } else if (anrm > rmax) {
  600. iscale = 1;
  601. sigma = rmax / anrm;
  602. }
  603. if (iscale == 1) {
  604. dlascl_(uplo, &c__0, &c__0, &c_b17, &sigma, n, n, &a[a_offset], lda,
  605. info);
  606. }
  607. /* Call DSYTRD to reduce symmetric matrix to tridiagonal form. */
  608. inde = 1;
  609. indtau = inde + *n;
  610. indwrk = indtau + *n;
  611. llwork = *lwork - indwrk + 1;
  612. dsytrd_(uplo, n, &a[a_offset], lda, &w[1], &work[inde], &work[indtau], &
  613. work[indwrk], &llwork, &iinfo);
  614. /* For eigenvalues only, call DSTERF. For eigenvectors, first call */
  615. /* DORGTR to generate the orthogonal matrix, then call DSTEQR. */
  616. if (! wantz) {
  617. dsterf_(n, &w[1], &work[inde], info);
  618. } else {
  619. dorgtr_(uplo, n, &a[a_offset], lda, &work[indtau], &work[indwrk], &
  620. llwork, &iinfo);
  621. dsteqr_(jobz, n, &w[1], &work[inde], &a[a_offset], lda, &work[indtau],
  622. info);
  623. }
  624. /* If matrix was scaled, then rescale eigenvalues appropriately. */
  625. if (iscale == 1) {
  626. if (*info == 0) {
  627. imax = *n;
  628. } else {
  629. imax = *info - 1;
  630. }
  631. d__1 = 1. / sigma;
  632. dscal_(&imax, &d__1, &w[1], &c__1);
  633. }
  634. /* Set WORK(1) to optimal workspace size. */
  635. work[1] = (doublereal) lwkopt;
  636. return 0;
  637. /* End of DSYEV */
  638. } /* dsyev_ */