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.

sgeev.c 28 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901
  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. #if defined(_WIN64)
  13. typedef long long BLASLONG;
  14. typedef unsigned long long BLASULONG;
  15. #else
  16. typedef long BLASLONG;
  17. typedef unsigned long BLASULONG;
  18. #endif
  19. #ifdef LAPACK_ILP64
  20. typedef BLASLONG blasint;
  21. #if defined(_WIN64)
  22. #define blasabs(x) llabs(x)
  23. #else
  24. #define blasabs(x) labs(x)
  25. #endif
  26. #else
  27. typedef int blasint;
  28. #define blasabs(x) abs(x)
  29. #endif
  30. typedef blasint integer;
  31. typedef unsigned int uinteger;
  32. typedef char *address;
  33. typedef short int shortint;
  34. typedef float real;
  35. typedef double doublereal;
  36. typedef struct { real r, i; } complex;
  37. typedef struct { doublereal r, i; } doublecomplex;
  38. #ifdef _MSC_VER
  39. static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
  40. static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
  41. static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
  42. static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
  43. #else
  44. static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
  45. static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
  46. static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
  47. static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
  48. #endif
  49. #define pCf(z) (*_pCf(z))
  50. #define pCd(z) (*_pCd(z))
  51. typedef blasint logical;
  52. typedef char logical1;
  53. typedef char integer1;
  54. #define TRUE_ (1)
  55. #define FALSE_ (0)
  56. /* Extern is for use with -E */
  57. #ifndef Extern
  58. #define Extern extern
  59. #endif
  60. /* I/O stuff */
  61. typedef int flag;
  62. typedef int ftnlen;
  63. typedef int ftnint;
  64. /*external read, write*/
  65. typedef struct
  66. { flag cierr;
  67. ftnint ciunit;
  68. flag ciend;
  69. char *cifmt;
  70. ftnint cirec;
  71. } cilist;
  72. /*internal read, write*/
  73. typedef struct
  74. { flag icierr;
  75. char *iciunit;
  76. flag iciend;
  77. char *icifmt;
  78. ftnint icirlen;
  79. ftnint icirnum;
  80. } icilist;
  81. /*open*/
  82. typedef struct
  83. { flag oerr;
  84. ftnint ounit;
  85. char *ofnm;
  86. ftnlen ofnmlen;
  87. char *osta;
  88. char *oacc;
  89. char *ofm;
  90. ftnint orl;
  91. char *oblnk;
  92. } olist;
  93. /*close*/
  94. typedef struct
  95. { flag cerr;
  96. ftnint cunit;
  97. char *csta;
  98. } cllist;
  99. /*rewind, backspace, endfile*/
  100. typedef struct
  101. { flag aerr;
  102. ftnint aunit;
  103. } alist;
  104. /* inquire */
  105. typedef struct
  106. { flag inerr;
  107. ftnint inunit;
  108. char *infile;
  109. ftnlen infilen;
  110. ftnint *inex; /*parameters in standard's order*/
  111. ftnint *inopen;
  112. ftnint *innum;
  113. ftnint *innamed;
  114. char *inname;
  115. ftnlen innamlen;
  116. char *inacc;
  117. ftnlen inacclen;
  118. char *inseq;
  119. ftnlen inseqlen;
  120. char *indir;
  121. ftnlen indirlen;
  122. char *infmt;
  123. ftnlen infmtlen;
  124. char *inform;
  125. ftnint informlen;
  126. char *inunf;
  127. ftnlen inunflen;
  128. ftnint *inrecl;
  129. ftnint *innrec;
  130. char *inblank;
  131. ftnlen inblanklen;
  132. } inlist;
  133. #define VOID void
  134. union Multitype { /* for multiple entry points */
  135. integer1 g;
  136. shortint h;
  137. integer i;
  138. /* longint j; */
  139. real r;
  140. doublereal d;
  141. complex c;
  142. doublecomplex z;
  143. };
  144. typedef union Multitype Multitype;
  145. struct Vardesc { /* for Namelist */
  146. char *name;
  147. char *addr;
  148. ftnlen *dims;
  149. int type;
  150. };
  151. typedef struct Vardesc Vardesc;
  152. struct Namelist {
  153. char *name;
  154. Vardesc **vars;
  155. int nvars;
  156. };
  157. typedef struct Namelist Namelist;
  158. #define abs(x) ((x) >= 0 ? (x) : -(x))
  159. #define dabs(x) (fabs(x))
  160. #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
  161. #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
  162. #define dmin(a,b) (f2cmin(a,b))
  163. #define dmax(a,b) (f2cmax(a,b))
  164. #define bit_test(a,b) ((a) >> (b) & 1)
  165. #define bit_clear(a,b) ((a) & ~((uinteger)1 << (b)))
  166. #define bit_set(a,b) ((a) | ((uinteger)1 << (b)))
  167. #define abort_() { sig_die("Fortran abort routine called", 1); }
  168. #define c_abs(z) (cabsf(Cf(z)))
  169. #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
  170. #ifdef _MSC_VER
  171. #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]);}
  172. #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]);}
  173. #else
  174. #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
  175. #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
  176. #endif
  177. #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
  178. #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
  179. #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
  180. //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
  181. #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
  182. #define d_abs(x) (fabs(*(x)))
  183. #define d_acos(x) (acos(*(x)))
  184. #define d_asin(x) (asin(*(x)))
  185. #define d_atan(x) (atan(*(x)))
  186. #define d_atn2(x, y) (atan2(*(x),*(y)))
  187. #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
  188. #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
  189. #define d_cos(x) (cos(*(x)))
  190. #define d_cosh(x) (cosh(*(x)))
  191. #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
  192. #define d_exp(x) (exp(*(x)))
  193. #define d_imag(z) (cimag(Cd(z)))
  194. #define r_imag(z) (cimagf(Cf(z)))
  195. #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
  196. #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
  197. #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
  198. #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
  199. #define d_log(x) (log(*(x)))
  200. #define d_mod(x, y) (fmod(*(x), *(y)))
  201. #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
  202. #define d_nint(x) u_nint(*(x))
  203. #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
  204. #define d_sign(a,b) u_sign(*(a),*(b))
  205. #define r_sign(a,b) u_sign(*(a),*(b))
  206. #define d_sin(x) (sin(*(x)))
  207. #define d_sinh(x) (sinh(*(x)))
  208. #define d_sqrt(x) (sqrt(*(x)))
  209. #define d_tan(x) (tan(*(x)))
  210. #define d_tanh(x) (tanh(*(x)))
  211. #define i_abs(x) abs(*(x))
  212. #define i_dnnt(x) ((integer)u_nint(*(x)))
  213. #define i_len(s, n) (n)
  214. #define i_nint(x) ((integer)u_nint(*(x)))
  215. #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
  216. #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++ = ' '; }
  217. #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
  218. #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]; }
  219. #define sig_die(s, kill) { exit(1); }
  220. #define s_stop(s, n) {exit(0);}
  221. #define z_abs(z) (cabs(Cd(z)))
  222. #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
  223. #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
  224. #define myexit_() break;
  225. #define mycycle() continue;
  226. #define myceiling(w) {ceil(w)}
  227. #define myhuge(w) {HUGE_VAL}
  228. //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
  229. #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
  230. /* -- translated by f2c (version 20000121).
  231. You must link the resulting object file with the libraries:
  232. -lf2c -lm (in that order)
  233. */
  234. /* Table of constant values */
  235. static integer c__1 = 1;
  236. static integer c__0 = 0;
  237. static integer c_n1 = -1;
  238. /* > \brief <b> SGEEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matr
  239. ices</b> */
  240. /* =========== DOCUMENTATION =========== */
  241. /* Online html documentation available at */
  242. /* http://www.netlib.org/lapack/explore-html/ */
  243. /* > \htmlonly */
  244. /* > Download SGEEV + dependencies */
  245. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgeev.f
  246. "> */
  247. /* > [TGZ]</a> */
  248. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgeev.f
  249. "> */
  250. /* > [ZIP]</a> */
  251. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgeev.f
  252. "> */
  253. /* > [TXT]</a> */
  254. /* > \endhtmlonly */
  255. /* Definition: */
  256. /* =========== */
  257. /* SUBROUTINE SGEEV( JOBVL, JOBVR, N, A, LDA, WR, WI, VL, LDVL, VR, */
  258. /* LDVR, WORK, LWORK, INFO ) */
  259. /* CHARACTER JOBVL, JOBVR */
  260. /* INTEGER INFO, LDA, LDVL, LDVR, LWORK, N */
  261. /* REAL A( LDA, * ), VL( LDVL, * ), VR( LDVR, * ), */
  262. /* $ WI( * ), WORK( * ), WR( * ) */
  263. /* > \par Purpose: */
  264. /* ============= */
  265. /* > */
  266. /* > \verbatim */
  267. /* > */
  268. /* > SGEEV computes for an N-by-N real nonsymmetric matrix A, the */
  269. /* > eigenvalues and, optionally, the left and/or right eigenvectors. */
  270. /* > */
  271. /* > The right eigenvector v(j) of A satisfies */
  272. /* > A * v(j) = lambda(j) * v(j) */
  273. /* > where lambda(j) is its eigenvalue. */
  274. /* > The left eigenvector u(j) of A satisfies */
  275. /* > u(j)**H * A = lambda(j) * u(j)**H */
  276. /* > where u(j)**H denotes the conjugate-transpose of u(j). */
  277. /* > */
  278. /* > The computed eigenvectors are normalized to have Euclidean norm */
  279. /* > equal to 1 and largest component real. */
  280. /* > \endverbatim */
  281. /* Arguments: */
  282. /* ========== */
  283. /* > \param[in] JOBVL */
  284. /* > \verbatim */
  285. /* > JOBVL is CHARACTER*1 */
  286. /* > = 'N': left eigenvectors of A are not computed; */
  287. /* > = 'V': left eigenvectors of A are computed. */
  288. /* > \endverbatim */
  289. /* > */
  290. /* > \param[in] JOBVR */
  291. /* > \verbatim */
  292. /* > JOBVR is CHARACTER*1 */
  293. /* > = 'N': right eigenvectors of A are not computed; */
  294. /* > = 'V': right eigenvectors of A are computed. */
  295. /* > \endverbatim */
  296. /* > */
  297. /* > \param[in] N */
  298. /* > \verbatim */
  299. /* > N is INTEGER */
  300. /* > The order of the matrix A. N >= 0. */
  301. /* > \endverbatim */
  302. /* > */
  303. /* > \param[in,out] A */
  304. /* > \verbatim */
  305. /* > A is REAL array, dimension (LDA,N) */
  306. /* > On entry, the N-by-N matrix A. */
  307. /* > On exit, A has been overwritten. */
  308. /* > \endverbatim */
  309. /* > */
  310. /* > \param[in] LDA */
  311. /* > \verbatim */
  312. /* > LDA is INTEGER */
  313. /* > The leading dimension of the array A. LDA >= f2cmax(1,N). */
  314. /* > \endverbatim */
  315. /* > */
  316. /* > \param[out] WR */
  317. /* > \verbatim */
  318. /* > WR is REAL array, dimension (N) */
  319. /* > \endverbatim */
  320. /* > */
  321. /* > \param[out] WI */
  322. /* > \verbatim */
  323. /* > WI is REAL array, dimension (N) */
  324. /* > WR and WI contain the real and imaginary parts, */
  325. /* > respectively, of the computed eigenvalues. Complex */
  326. /* > conjugate pairs of eigenvalues appear consecutively */
  327. /* > with the eigenvalue having the positive imaginary part */
  328. /* > first. */
  329. /* > \endverbatim */
  330. /* > */
  331. /* > \param[out] VL */
  332. /* > \verbatim */
  333. /* > VL is REAL array, dimension (LDVL,N) */
  334. /* > If JOBVL = 'V', the left eigenvectors u(j) are stored one */
  335. /* > after another in the columns of VL, in the same order */
  336. /* > as their eigenvalues. */
  337. /* > If JOBVL = 'N', VL is not referenced. */
  338. /* > If the j-th eigenvalue is real, then u(j) = VL(:,j), */
  339. /* > the j-th column of VL. */
  340. /* > If the j-th and (j+1)-st eigenvalues form a complex */
  341. /* > conjugate pair, then u(j) = VL(:,j) + i*VL(:,j+1) and */
  342. /* > u(j+1) = VL(:,j) - i*VL(:,j+1). */
  343. /* > \endverbatim */
  344. /* > */
  345. /* > \param[in] LDVL */
  346. /* > \verbatim */
  347. /* > LDVL is INTEGER */
  348. /* > The leading dimension of the array VL. LDVL >= 1; if */
  349. /* > JOBVL = 'V', LDVL >= N. */
  350. /* > \endverbatim */
  351. /* > */
  352. /* > \param[out] VR */
  353. /* > \verbatim */
  354. /* > VR is REAL array, dimension (LDVR,N) */
  355. /* > If JOBVR = 'V', the right eigenvectors v(j) are stored one */
  356. /* > after another in the columns of VR, in the same order */
  357. /* > as their eigenvalues. */
  358. /* > If JOBVR = 'N', VR is not referenced. */
  359. /* > If the j-th eigenvalue is real, then v(j) = VR(:,j), */
  360. /* > the j-th column of VR. */
  361. /* > If the j-th and (j+1)-st eigenvalues form a complex */
  362. /* > conjugate pair, then v(j) = VR(:,j) + i*VR(:,j+1) and */
  363. /* > v(j+1) = VR(:,j) - i*VR(:,j+1). */
  364. /* > \endverbatim */
  365. /* > */
  366. /* > \param[in] LDVR */
  367. /* > \verbatim */
  368. /* > LDVR is INTEGER */
  369. /* > The leading dimension of the array VR. LDVR >= 1; if */
  370. /* > JOBVR = 'V', LDVR >= N. */
  371. /* > \endverbatim */
  372. /* > */
  373. /* > \param[out] WORK */
  374. /* > \verbatim */
  375. /* > WORK is REAL array, dimension (MAX(1,LWORK)) */
  376. /* > On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
  377. /* > \endverbatim */
  378. /* > */
  379. /* > \param[in] LWORK */
  380. /* > \verbatim */
  381. /* > LWORK is INTEGER */
  382. /* > The dimension of the array WORK. LWORK >= f2cmax(1,3*N), and */
  383. /* > if JOBVL = 'V' or JOBVR = 'V', LWORK >= 4*N. For good */
  384. /* > performance, LWORK must generally be larger. */
  385. /* > */
  386. /* > If LWORK = -1, then a workspace query is assumed; the routine */
  387. /* > only calculates the optimal size of the WORK array, returns */
  388. /* > this value as the first entry of the WORK array, and no error */
  389. /* > message related to LWORK is issued by XERBLA. */
  390. /* > \endverbatim */
  391. /* > */
  392. /* > \param[out] INFO */
  393. /* > \verbatim */
  394. /* > INFO is INTEGER */
  395. /* > = 0: successful exit */
  396. /* > < 0: if INFO = -i, the i-th argument had an illegal value. */
  397. /* > > 0: if INFO = i, the QR algorithm failed to compute all the */
  398. /* > eigenvalues, and no eigenvectors have been computed; */
  399. /* > elements i+1:N of WR and WI contain eigenvalues which */
  400. /* > have converged. */
  401. /* > \endverbatim */
  402. /* Authors: */
  403. /* ======== */
  404. /* > \author Univ. of Tennessee */
  405. /* > \author Univ. of California Berkeley */
  406. /* > \author Univ. of Colorado Denver */
  407. /* > \author NAG Ltd. */
  408. /* > \date June 2016 */
  409. /* @generated from dgeev.f, fortran d -> s, Tue Apr 19 01:47:44 2016 */
  410. /* > \ingroup realGEeigen */
  411. /* ===================================================================== */
  412. /* Subroutine */ void sgeev_(char *jobvl, char *jobvr, integer *n, real *a,
  413. integer *lda, real *wr, real *wi, real *vl, integer *ldvl, real *vr,
  414. integer *ldvr, real *work, integer *lwork, integer *info)
  415. {
  416. /* System generated locals */
  417. integer a_dim1, a_offset, vl_dim1, vl_offset, vr_dim1, vr_offset, i__1,
  418. i__2, i__3;
  419. real r__1, r__2;
  420. /* Local variables */
  421. integer ibal;
  422. char side[1];
  423. real anrm;
  424. integer ierr, itau, iwrk, nout;
  425. extern /* Subroutine */ void srot_(integer *, real *, integer *, real *,
  426. integer *, real *, real *);
  427. extern real snrm2_(integer *, real *, integer *);
  428. integer i__, k;
  429. real r__;
  430. extern logical lsame_(char *, char *);
  431. extern /* Subroutine */ void sscal_(integer *, real *, real *, integer *);
  432. extern real slapy2_(real *, real *);
  433. real cs;
  434. extern /* Subroutine */ void slabad_(real *, real *);
  435. logical scalea;
  436. real cscale;
  437. extern /* Subroutine */ void sgebak_(char *, char *, integer *, integer *,
  438. integer *, real *, integer *, real *, integer *, integer *), sgebal_(char *, integer *, real *, integer *,
  439. integer *, integer *, real *, integer *);
  440. real sn;
  441. extern real slamch_(char *), slange_(char *, integer *, integer *,
  442. real *, integer *, real *);
  443. extern /* Subroutine */ void sgehrd_(integer *, integer *, integer *, real
  444. *, integer *, real *, real *, integer *, integer *);
  445. extern int xerbla_(char *, integer *, ftnlen);
  446. extern integer ilaenv_(integer *, char *, char *, integer *, integer *,
  447. integer *, integer *, ftnlen, ftnlen);
  448. logical select[1];
  449. real bignum;
  450. extern /* Subroutine */ void slascl_(char *, integer *, integer *, real *,
  451. real *, integer *, integer *, real *, integer *, integer *);
  452. extern integer isamax_(integer *, real *, integer *);
  453. extern /* Subroutine */ void slacpy_(char *, integer *, integer *, real *,
  454. integer *, real *, integer *), slartg_(real *, real *,
  455. real *, real *, real *), sorghr_(integer *, integer *, integer *,
  456. real *, integer *, real *, real *, integer *, integer *), shseqr_(
  457. char *, char *, integer *, integer *, integer *, real *, integer *
  458. , real *, real *, real *, integer *, real *, integer *, integer *);
  459. integer minwrk, maxwrk;
  460. logical wantvl;
  461. real smlnum;
  462. integer hswork;
  463. logical lquery, wantvr;
  464. extern /* Subroutine */ void strevc3_(char *, char *, logical *, integer *,
  465. real *, integer *, real *, integer *, real *, integer *, integer
  466. *, integer *, real *, integer *, integer *);
  467. integer ihi;
  468. real scl;
  469. integer ilo;
  470. real dum[1], eps;
  471. integer lwork_trevc__;
  472. /* -- LAPACK driver routine (version 3.7.0) -- */
  473. /* -- LAPACK is a software package provided by Univ. of Tennessee, -- */
  474. /* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
  475. /* June 2016 */
  476. /* ===================================================================== */
  477. /* Test the input arguments */
  478. /* Parameter adjustments */
  479. a_dim1 = *lda;
  480. a_offset = 1 + a_dim1 * 1;
  481. a -= a_offset;
  482. --wr;
  483. --wi;
  484. vl_dim1 = *ldvl;
  485. vl_offset = 1 + vl_dim1 * 1;
  486. vl -= vl_offset;
  487. vr_dim1 = *ldvr;
  488. vr_offset = 1 + vr_dim1 * 1;
  489. vr -= vr_offset;
  490. --work;
  491. /* Function Body */
  492. *info = 0;
  493. lquery = *lwork == -1;
  494. wantvl = lsame_(jobvl, "V");
  495. wantvr = lsame_(jobvr, "V");
  496. if (! wantvl && ! lsame_(jobvl, "N")) {
  497. *info = -1;
  498. } else if (! wantvr && ! lsame_(jobvr, "N")) {
  499. *info = -2;
  500. } else if (*n < 0) {
  501. *info = -3;
  502. } else if (*lda < f2cmax(1,*n)) {
  503. *info = -5;
  504. } else if (*ldvl < 1 || wantvl && *ldvl < *n) {
  505. *info = -9;
  506. } else if (*ldvr < 1 || wantvr && *ldvr < *n) {
  507. *info = -11;
  508. }
  509. /* Compute workspace */
  510. /* (Note: Comments in the code beginning "Workspace:" describe the */
  511. /* minimal amount of workspace needed at that point in the code, */
  512. /* as well as the preferred amount for good performance. */
  513. /* NB refers to the optimal block size for the immediately */
  514. /* following subroutine, as returned by ILAENV. */
  515. /* HSWORK refers to the workspace preferred by SHSEQR, as */
  516. /* calculated below. HSWORK is computed assuming ILO=1 and IHI=N, */
  517. /* the worst case.) */
  518. if (*info == 0) {
  519. if (*n == 0) {
  520. minwrk = 1;
  521. maxwrk = 1;
  522. } else {
  523. maxwrk = (*n << 1) + *n * ilaenv_(&c__1, "SGEHRD", " ", n, &c__1,
  524. n, &c__0, (ftnlen)6, (ftnlen)1);
  525. if (wantvl) {
  526. minwrk = *n << 2;
  527. /* Computing MAX */
  528. i__1 = maxwrk, i__2 = (*n << 1) + (*n - 1) * ilaenv_(&c__1,
  529. "SORGHR", " ", n, &c__1, n, &c_n1, (ftnlen)6, (ftnlen)
  530. 1);
  531. maxwrk = f2cmax(i__1,i__2);
  532. shseqr_("S", "V", n, &c__1, n, &a[a_offset], lda, &wr[1], &wi[
  533. 1], &vl[vl_offset], ldvl, &work[1], &c_n1, info);
  534. hswork = (integer) work[1];
  535. /* Computing MAX */
  536. i__1 = maxwrk, i__2 = *n + 1, i__1 = f2cmax(i__1,i__2), i__2 = *
  537. n + hswork;
  538. maxwrk = f2cmax(i__1,i__2);
  539. strevc3_("L", "B", select, n, &a[a_offset], lda, &vl[
  540. vl_offset], ldvl, &vr[vr_offset], ldvr, n, &nout, &
  541. work[1], &c_n1, &ierr);
  542. lwork_trevc__ = (integer) work[1];
  543. /* Computing MAX */
  544. i__1 = maxwrk, i__2 = *n + lwork_trevc__;
  545. maxwrk = f2cmax(i__1,i__2);
  546. /* Computing MAX */
  547. i__1 = maxwrk, i__2 = *n << 2;
  548. maxwrk = f2cmax(i__1,i__2);
  549. } else if (wantvr) {
  550. minwrk = *n << 2;
  551. /* Computing MAX */
  552. i__1 = maxwrk, i__2 = (*n << 1) + (*n - 1) * ilaenv_(&c__1,
  553. "SORGHR", " ", n, &c__1, n, &c_n1, (ftnlen)6, (ftnlen)
  554. 1);
  555. maxwrk = f2cmax(i__1,i__2);
  556. shseqr_("S", "V", n, &c__1, n, &a[a_offset], lda, &wr[1], &wi[
  557. 1], &vr[vr_offset], ldvr, &work[1], &c_n1, info);
  558. hswork = (integer) work[1];
  559. /* Computing MAX */
  560. i__1 = maxwrk, i__2 = *n + 1, i__1 = f2cmax(i__1,i__2), i__2 = *
  561. n + hswork;
  562. maxwrk = f2cmax(i__1,i__2);
  563. strevc3_("R", "B", select, n, &a[a_offset], lda, &vl[
  564. vl_offset], ldvl, &vr[vr_offset], ldvr, n, &nout, &
  565. work[1], &c_n1, &ierr);
  566. lwork_trevc__ = (integer) work[1];
  567. /* Computing MAX */
  568. i__1 = maxwrk, i__2 = *n + lwork_trevc__;
  569. maxwrk = f2cmax(i__1,i__2);
  570. /* Computing MAX */
  571. i__1 = maxwrk, i__2 = *n << 2;
  572. maxwrk = f2cmax(i__1,i__2);
  573. } else {
  574. minwrk = *n * 3;
  575. shseqr_("E", "N", n, &c__1, n, &a[a_offset], lda, &wr[1], &wi[
  576. 1], &vr[vr_offset], ldvr, &work[1], &c_n1, info);
  577. hswork = (integer) work[1];
  578. /* Computing MAX */
  579. i__1 = maxwrk, i__2 = *n + 1, i__1 = f2cmax(i__1,i__2), i__2 = *
  580. n + hswork;
  581. maxwrk = f2cmax(i__1,i__2);
  582. }
  583. maxwrk = f2cmax(maxwrk,minwrk);
  584. }
  585. work[1] = (real) maxwrk;
  586. if (*lwork < minwrk && ! lquery) {
  587. *info = -13;
  588. }
  589. }
  590. if (*info != 0) {
  591. i__1 = -(*info);
  592. xerbla_("SGEEV ", &i__1, (ftnlen)5);
  593. return;
  594. } else if (lquery) {
  595. return;
  596. }
  597. /* Quick return if possible */
  598. if (*n == 0) {
  599. return;
  600. }
  601. /* Get machine constants */
  602. eps = slamch_("P");
  603. smlnum = slamch_("S");
  604. bignum = 1.f / smlnum;
  605. slabad_(&smlnum, &bignum);
  606. smlnum = sqrt(smlnum) / eps;
  607. bignum = 1.f / smlnum;
  608. /* Scale A if f2cmax element outside range [SMLNUM,BIGNUM] */
  609. anrm = slange_("M", n, n, &a[a_offset], lda, dum);
  610. scalea = FALSE_;
  611. if (anrm > 0.f && anrm < smlnum) {
  612. scalea = TRUE_;
  613. cscale = smlnum;
  614. } else if (anrm > bignum) {
  615. scalea = TRUE_;
  616. cscale = bignum;
  617. }
  618. if (scalea) {
  619. slascl_("G", &c__0, &c__0, &anrm, &cscale, n, n, &a[a_offset], lda, &
  620. ierr);
  621. }
  622. /* Balance the matrix */
  623. /* (Workspace: need N) */
  624. ibal = 1;
  625. sgebal_("B", n, &a[a_offset], lda, &ilo, &ihi, &work[ibal], &ierr);
  626. /* Reduce to upper Hessenberg form */
  627. /* (Workspace: need 3*N, prefer 2*N+N*NB) */
  628. itau = ibal + *n;
  629. iwrk = itau + *n;
  630. i__1 = *lwork - iwrk + 1;
  631. sgehrd_(n, &ilo, &ihi, &a[a_offset], lda, &work[itau], &work[iwrk], &i__1,
  632. &ierr);
  633. if (wantvl) {
  634. /* Want left eigenvectors */
  635. /* Copy Householder vectors to VL */
  636. *(unsigned char *)side = 'L';
  637. slacpy_("L", n, n, &a[a_offset], lda, &vl[vl_offset], ldvl)
  638. ;
  639. /* Generate orthogonal matrix in VL */
  640. /* (Workspace: need 3*N-1, prefer 2*N+(N-1)*NB) */
  641. i__1 = *lwork - iwrk + 1;
  642. sorghr_(n, &ilo, &ihi, &vl[vl_offset], ldvl, &work[itau], &work[iwrk],
  643. &i__1, &ierr);
  644. /* Perform QR iteration, accumulating Schur vectors in VL */
  645. /* (Workspace: need N+1, prefer N+HSWORK (see comments) ) */
  646. iwrk = itau;
  647. i__1 = *lwork - iwrk + 1;
  648. shseqr_("S", "V", n, &ilo, &ihi, &a[a_offset], lda, &wr[1], &wi[1], &
  649. vl[vl_offset], ldvl, &work[iwrk], &i__1, info);
  650. if (wantvr) {
  651. /* Want left and right eigenvectors */
  652. /* Copy Schur vectors to VR */
  653. *(unsigned char *)side = 'B';
  654. slacpy_("F", n, n, &vl[vl_offset], ldvl, &vr[vr_offset], ldvr);
  655. }
  656. } else if (wantvr) {
  657. /* Want right eigenvectors */
  658. /* Copy Householder vectors to VR */
  659. *(unsigned char *)side = 'R';
  660. slacpy_("L", n, n, &a[a_offset], lda, &vr[vr_offset], ldvr)
  661. ;
  662. /* Generate orthogonal matrix in VR */
  663. /* (Workspace: need 3*N-1, prefer 2*N+(N-1)*NB) */
  664. i__1 = *lwork - iwrk + 1;
  665. sorghr_(n, &ilo, &ihi, &vr[vr_offset], ldvr, &work[itau], &work[iwrk],
  666. &i__1, &ierr);
  667. /* Perform QR iteration, accumulating Schur vectors in VR */
  668. /* (Workspace: need N+1, prefer N+HSWORK (see comments) ) */
  669. iwrk = itau;
  670. i__1 = *lwork - iwrk + 1;
  671. shseqr_("S", "V", n, &ilo, &ihi, &a[a_offset], lda, &wr[1], &wi[1], &
  672. vr[vr_offset], ldvr, &work[iwrk], &i__1, info);
  673. } else {
  674. /* Compute eigenvalues only */
  675. /* (Workspace: need N+1, prefer N+HSWORK (see comments) ) */
  676. iwrk = itau;
  677. i__1 = *lwork - iwrk + 1;
  678. shseqr_("E", "N", n, &ilo, &ihi, &a[a_offset], lda, &wr[1], &wi[1], &
  679. vr[vr_offset], ldvr, &work[iwrk], &i__1, info);
  680. }
  681. /* If INFO .NE. 0 from SHSEQR, then quit */
  682. if (*info != 0) {
  683. goto L50;
  684. }
  685. if (wantvl || wantvr) {
  686. /* Compute left and/or right eigenvectors */
  687. /* (Workspace: need 4*N, prefer N + N + 2*N*NB) */
  688. i__1 = *lwork - iwrk + 1;
  689. strevc3_(side, "B", select, n, &a[a_offset], lda, &vl[vl_offset],
  690. ldvl, &vr[vr_offset], ldvr, n, &nout, &work[iwrk], &i__1, &
  691. ierr);
  692. }
  693. if (wantvl) {
  694. /* Undo balancing of left eigenvectors */
  695. /* (Workspace: need N) */
  696. sgebak_("B", "L", n, &ilo, &ihi, &work[ibal], n, &vl[vl_offset], ldvl,
  697. &ierr);
  698. /* Normalize left eigenvectors and make largest component real */
  699. i__1 = *n;
  700. for (i__ = 1; i__ <= i__1; ++i__) {
  701. if (wi[i__] == 0.f) {
  702. scl = 1.f / snrm2_(n, &vl[i__ * vl_dim1 + 1], &c__1);
  703. sscal_(n, &scl, &vl[i__ * vl_dim1 + 1], &c__1);
  704. } else if (wi[i__] > 0.f) {
  705. r__1 = snrm2_(n, &vl[i__ * vl_dim1 + 1], &c__1);
  706. r__2 = snrm2_(n, &vl[(i__ + 1) * vl_dim1 + 1], &c__1);
  707. scl = 1.f / slapy2_(&r__1, &r__2);
  708. sscal_(n, &scl, &vl[i__ * vl_dim1 + 1], &c__1);
  709. sscal_(n, &scl, &vl[(i__ + 1) * vl_dim1 + 1], &c__1);
  710. i__2 = *n;
  711. for (k = 1; k <= i__2; ++k) {
  712. /* Computing 2nd power */
  713. r__1 = vl[k + i__ * vl_dim1];
  714. /* Computing 2nd power */
  715. r__2 = vl[k + (i__ + 1) * vl_dim1];
  716. work[iwrk + k - 1] = r__1 * r__1 + r__2 * r__2;
  717. /* L10: */
  718. }
  719. k = isamax_(n, &work[iwrk], &c__1);
  720. slartg_(&vl[k + i__ * vl_dim1], &vl[k + (i__ + 1) * vl_dim1],
  721. &cs, &sn, &r__);
  722. srot_(n, &vl[i__ * vl_dim1 + 1], &c__1, &vl[(i__ + 1) *
  723. vl_dim1 + 1], &c__1, &cs, &sn);
  724. vl[k + (i__ + 1) * vl_dim1] = 0.f;
  725. }
  726. /* L20: */
  727. }
  728. }
  729. if (wantvr) {
  730. /* Undo balancing of right eigenvectors */
  731. /* (Workspace: need N) */
  732. sgebak_("B", "R", n, &ilo, &ihi, &work[ibal], n, &vr[vr_offset], ldvr,
  733. &ierr);
  734. /* Normalize right eigenvectors and make largest component real */
  735. i__1 = *n;
  736. for (i__ = 1; i__ <= i__1; ++i__) {
  737. if (wi[i__] == 0.f) {
  738. scl = 1.f / snrm2_(n, &vr[i__ * vr_dim1 + 1], &c__1);
  739. sscal_(n, &scl, &vr[i__ * vr_dim1 + 1], &c__1);
  740. } else if (wi[i__] > 0.f) {
  741. r__1 = snrm2_(n, &vr[i__ * vr_dim1 + 1], &c__1);
  742. r__2 = snrm2_(n, &vr[(i__ + 1) * vr_dim1 + 1], &c__1);
  743. scl = 1.f / slapy2_(&r__1, &r__2);
  744. sscal_(n, &scl, &vr[i__ * vr_dim1 + 1], &c__1);
  745. sscal_(n, &scl, &vr[(i__ + 1) * vr_dim1 + 1], &c__1);
  746. i__2 = *n;
  747. for (k = 1; k <= i__2; ++k) {
  748. /* Computing 2nd power */
  749. r__1 = vr[k + i__ * vr_dim1];
  750. /* Computing 2nd power */
  751. r__2 = vr[k + (i__ + 1) * vr_dim1];
  752. work[iwrk + k - 1] = r__1 * r__1 + r__2 * r__2;
  753. /* L30: */
  754. }
  755. k = isamax_(n, &work[iwrk], &c__1);
  756. slartg_(&vr[k + i__ * vr_dim1], &vr[k + (i__ + 1) * vr_dim1],
  757. &cs, &sn, &r__);
  758. srot_(n, &vr[i__ * vr_dim1 + 1], &c__1, &vr[(i__ + 1) *
  759. vr_dim1 + 1], &c__1, &cs, &sn);
  760. vr[k + (i__ + 1) * vr_dim1] = 0.f;
  761. }
  762. /* L40: */
  763. }
  764. }
  765. /* Undo scaling if necessary */
  766. L50:
  767. if (scalea) {
  768. i__1 = *n - *info;
  769. /* Computing MAX */
  770. i__3 = *n - *info;
  771. i__2 = f2cmax(i__3,1);
  772. slascl_("G", &c__0, &c__0, &cscale, &anrm, &i__1, &c__1, &wr[*info +
  773. 1], &i__2, &ierr);
  774. i__1 = *n - *info;
  775. /* Computing MAX */
  776. i__3 = *n - *info;
  777. i__2 = f2cmax(i__3,1);
  778. slascl_("G", &c__0, &c__0, &cscale, &anrm, &i__1, &c__1, &wi[*info +
  779. 1], &i__2, &ierr);
  780. if (*info > 0) {
  781. i__1 = ilo - 1;
  782. slascl_("G", &c__0, &c__0, &cscale, &anrm, &i__1, &c__1, &wr[1],
  783. n, &ierr);
  784. i__1 = ilo - 1;
  785. slascl_("G", &c__0, &c__0, &cscale, &anrm, &i__1, &c__1, &wi[1],
  786. n, &ierr);
  787. }
  788. }
  789. work[1] = (real) maxwrk;
  790. return;
  791. /* End of SGEEV */
  792. } /* sgeev_ */