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.

sggev.c 28 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972
  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. static real c_b36 = 0.f;
  239. static real c_b37 = 1.f;
  240. /* > \brief <b> SGGEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matr
  241. ices</b> */
  242. /* =========== DOCUMENTATION =========== */
  243. /* Online html documentation available at */
  244. /* http://www.netlib.org/lapack/explore-html/ */
  245. /* > \htmlonly */
  246. /* > Download SGGEV + dependencies */
  247. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sggev.f
  248. "> */
  249. /* > [TGZ]</a> */
  250. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sggev.f
  251. "> */
  252. /* > [ZIP]</a> */
  253. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sggev.f
  254. "> */
  255. /* > [TXT]</a> */
  256. /* > \endhtmlonly */
  257. /* Definition: */
  258. /* =========== */
  259. /* SUBROUTINE SGGEV( JOBVL, JOBVR, N, A, LDA, B, LDB, ALPHAR, ALPHAI, */
  260. /* BETA, VL, LDVL, VR, LDVR, WORK, LWORK, INFO ) */
  261. /* CHARACTER JOBVL, JOBVR */
  262. /* INTEGER INFO, LDA, LDB, LDVL, LDVR, LWORK, N */
  263. /* REAL A( LDA, * ), ALPHAI( * ), ALPHAR( * ), */
  264. /* $ B( LDB, * ), BETA( * ), VL( LDVL, * ), */
  265. /* $ VR( LDVR, * ), WORK( * ) */
  266. /* > \par Purpose: */
  267. /* ============= */
  268. /* > */
  269. /* > \verbatim */
  270. /* > */
  271. /* > SGGEV computes for a pair of N-by-N real nonsymmetric matrices (A,B) */
  272. /* > the generalized eigenvalues, and optionally, the left and/or right */
  273. /* > generalized eigenvectors. */
  274. /* > */
  275. /* > A generalized eigenvalue for a pair of matrices (A,B) is a scalar */
  276. /* > lambda or a ratio alpha/beta = lambda, such that A - lambda*B is */
  277. /* > singular. It is usually represented as the pair (alpha,beta), as */
  278. /* > there is a reasonable interpretation for beta=0, and even for both */
  279. /* > being zero. */
  280. /* > */
  281. /* > The right eigenvector v(j) corresponding to the eigenvalue lambda(j) */
  282. /* > of (A,B) satisfies */
  283. /* > */
  284. /* > A * v(j) = lambda(j) * B * v(j). */
  285. /* > */
  286. /* > The left eigenvector u(j) corresponding to the eigenvalue lambda(j) */
  287. /* > of (A,B) satisfies */
  288. /* > */
  289. /* > u(j)**H * A = lambda(j) * u(j)**H * B . */
  290. /* > */
  291. /* > where u(j)**H is the conjugate-transpose of u(j). */
  292. /* > */
  293. /* > \endverbatim */
  294. /* Arguments: */
  295. /* ========== */
  296. /* > \param[in] JOBVL */
  297. /* > \verbatim */
  298. /* > JOBVL is CHARACTER*1 */
  299. /* > = 'N': do not compute the left generalized eigenvectors; */
  300. /* > = 'V': compute the left generalized eigenvectors. */
  301. /* > \endverbatim */
  302. /* > */
  303. /* > \param[in] JOBVR */
  304. /* > \verbatim */
  305. /* > JOBVR is CHARACTER*1 */
  306. /* > = 'N': do not compute the right generalized eigenvectors; */
  307. /* > = 'V': compute the right generalized eigenvectors. */
  308. /* > \endverbatim */
  309. /* > */
  310. /* > \param[in] N */
  311. /* > \verbatim */
  312. /* > N is INTEGER */
  313. /* > The order of the matrices A, B, VL, and VR. N >= 0. */
  314. /* > \endverbatim */
  315. /* > */
  316. /* > \param[in,out] A */
  317. /* > \verbatim */
  318. /* > A is REAL array, dimension (LDA, N) */
  319. /* > On entry, the matrix A in the pair (A,B). */
  320. /* > On exit, A has been overwritten. */
  321. /* > \endverbatim */
  322. /* > */
  323. /* > \param[in] LDA */
  324. /* > \verbatim */
  325. /* > LDA is INTEGER */
  326. /* > The leading dimension of A. LDA >= f2cmax(1,N). */
  327. /* > \endverbatim */
  328. /* > */
  329. /* > \param[in,out] B */
  330. /* > \verbatim */
  331. /* > B is REAL array, dimension (LDB, N) */
  332. /* > On entry, the matrix B in the pair (A,B). */
  333. /* > On exit, B has been overwritten. */
  334. /* > \endverbatim */
  335. /* > */
  336. /* > \param[in] LDB */
  337. /* > \verbatim */
  338. /* > LDB is INTEGER */
  339. /* > The leading dimension of B. LDB >= f2cmax(1,N). */
  340. /* > \endverbatim */
  341. /* > */
  342. /* > \param[out] ALPHAR */
  343. /* > \verbatim */
  344. /* > ALPHAR is REAL array, dimension (N) */
  345. /* > \endverbatim */
  346. /* > */
  347. /* > \param[out] ALPHAI */
  348. /* > \verbatim */
  349. /* > ALPHAI is REAL array, dimension (N) */
  350. /* > \endverbatim */
  351. /* > */
  352. /* > \param[out] BETA */
  353. /* > \verbatim */
  354. /* > BETA is REAL array, dimension (N) */
  355. /* > On exit, (ALPHAR(j) + ALPHAI(j)*i)/BETA(j), j=1,...,N, will */
  356. /* > be the generalized eigenvalues. If ALPHAI(j) is zero, then */
  357. /* > the j-th eigenvalue is real; if positive, then the j-th and */
  358. /* > (j+1)-st eigenvalues are a complex conjugate pair, with */
  359. /* > ALPHAI(j+1) negative. */
  360. /* > */
  361. /* > Note: the quotients ALPHAR(j)/BETA(j) and ALPHAI(j)/BETA(j) */
  362. /* > may easily over- or underflow, and BETA(j) may even be zero. */
  363. /* > Thus, the user should avoid naively computing the ratio */
  364. /* > alpha/beta. However, ALPHAR and ALPHAI will be always less */
  365. /* > than and usually comparable with norm(A) in magnitude, and */
  366. /* > BETA always less than and usually comparable with norm(B). */
  367. /* > \endverbatim */
  368. /* > */
  369. /* > \param[out] VL */
  370. /* > \verbatim */
  371. /* > VL is REAL array, dimension (LDVL,N) */
  372. /* > If JOBVL = 'V', the left eigenvectors u(j) are stored one */
  373. /* > after another in the columns of VL, in the same order as */
  374. /* > their eigenvalues. If the j-th eigenvalue is real, then */
  375. /* > u(j) = VL(:,j), the j-th column of VL. If the j-th and */
  376. /* > (j+1)-th eigenvalues form a complex conjugate pair, then */
  377. /* > u(j) = VL(:,j)+i*VL(:,j+1) and u(j+1) = VL(:,j)-i*VL(:,j+1). */
  378. /* > Each eigenvector is scaled so the largest component has */
  379. /* > abs(real part)+abs(imag. part)=1. */
  380. /* > Not referenced if JOBVL = 'N'. */
  381. /* > \endverbatim */
  382. /* > */
  383. /* > \param[in] LDVL */
  384. /* > \verbatim */
  385. /* > LDVL is INTEGER */
  386. /* > The leading dimension of the matrix VL. LDVL >= 1, and */
  387. /* > if JOBVL = 'V', LDVL >= N. */
  388. /* > \endverbatim */
  389. /* > */
  390. /* > \param[out] VR */
  391. /* > \verbatim */
  392. /* > VR is REAL array, dimension (LDVR,N) */
  393. /* > If JOBVR = 'V', the right eigenvectors v(j) are stored one */
  394. /* > after another in the columns of VR, in the same order as */
  395. /* > their eigenvalues. If the j-th eigenvalue is real, then */
  396. /* > v(j) = VR(:,j), the j-th column of VR. If the j-th and */
  397. /* > (j+1)-th eigenvalues form a complex conjugate pair, then */
  398. /* > v(j) = VR(:,j)+i*VR(:,j+1) and v(j+1) = VR(:,j)-i*VR(:,j+1). */
  399. /* > Each eigenvector is scaled so the largest component has */
  400. /* > abs(real part)+abs(imag. part)=1. */
  401. /* > Not referenced if JOBVR = 'N'. */
  402. /* > \endverbatim */
  403. /* > */
  404. /* > \param[in] LDVR */
  405. /* > \verbatim */
  406. /* > LDVR is INTEGER */
  407. /* > The leading dimension of the matrix VR. LDVR >= 1, and */
  408. /* > if JOBVR = 'V', LDVR >= N. */
  409. /* > \endverbatim */
  410. /* > */
  411. /* > \param[out] WORK */
  412. /* > \verbatim */
  413. /* > WORK is REAL array, dimension (MAX(1,LWORK)) */
  414. /* > On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
  415. /* > \endverbatim */
  416. /* > */
  417. /* > \param[in] LWORK */
  418. /* > \verbatim */
  419. /* > LWORK is INTEGER */
  420. /* > The dimension of the array WORK. LWORK >= f2cmax(1,8*N). */
  421. /* > For good performance, LWORK must generally be larger. */
  422. /* > */
  423. /* > If LWORK = -1, then a workspace query is assumed; the routine */
  424. /* > only calculates the optimal size of the WORK array, returns */
  425. /* > this value as the first entry of the WORK array, and no error */
  426. /* > message related to LWORK is issued by XERBLA. */
  427. /* > \endverbatim */
  428. /* > */
  429. /* > \param[out] INFO */
  430. /* > \verbatim */
  431. /* > INFO is INTEGER */
  432. /* > = 0: successful exit */
  433. /* > < 0: if INFO = -i, the i-th argument had an illegal value. */
  434. /* > = 1,...,N: */
  435. /* > The QZ iteration failed. No eigenvectors have been */
  436. /* > calculated, but ALPHAR(j), ALPHAI(j), and BETA(j) */
  437. /* > should be correct for j=INFO+1,...,N. */
  438. /* > > N: =N+1: other than QZ iteration failed in SHGEQZ. */
  439. /* > =N+2: error return from STGEVC. */
  440. /* > \endverbatim */
  441. /* Authors: */
  442. /* ======== */
  443. /* > \author Univ. of Tennessee */
  444. /* > \author Univ. of California Berkeley */
  445. /* > \author Univ. of Colorado Denver */
  446. /* > \author NAG Ltd. */
  447. /* > \date April 2012 */
  448. /* > \ingroup realGEeigen */
  449. /* ===================================================================== */
  450. /* Subroutine */ void sggev_(char *jobvl, char *jobvr, integer *n, real *a,
  451. integer *lda, real *b, integer *ldb, real *alphar, real *alphai, real
  452. *beta, real *vl, integer *ldvl, real *vr, integer *ldvr, real *work,
  453. integer *lwork, integer *info)
  454. {
  455. /* System generated locals */
  456. integer a_dim1, a_offset, b_dim1, b_offset, vl_dim1, vl_offset, vr_dim1,
  457. vr_offset, i__1, i__2;
  458. real r__1, r__2, r__3, r__4;
  459. /* Local variables */
  460. real anrm, bnrm;
  461. integer ierr, itau;
  462. real temp;
  463. logical ilvl, ilvr;
  464. integer iwrk;
  465. extern logical lsame_(char *, char *);
  466. integer ileft, icols, irows, jc;
  467. extern /* Subroutine */ void slabad_(real *, real *);
  468. integer in, jr;
  469. extern /* Subroutine */ void sggbak_(char *, char *, integer *, integer *,
  470. integer *, real *, real *, integer *, real *, integer *, integer *
  471. ), sggbal_(char *, integer *, real *, integer *,
  472. real *, integer *, integer *, integer *, real *, real *, real *,
  473. integer *);
  474. logical ilascl, ilbscl;
  475. extern real slamch_(char *), slange_(char *, integer *, integer *,
  476. real *, integer *, real *);
  477. extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
  478. extern void sgghrd_(
  479. char *, char *, integer *, integer *, integer *, real *, integer *
  480. , real *, integer *, real *, integer *, real *, integer *,
  481. integer *);
  482. logical ldumma[1];
  483. char chtemp[1];
  484. real bignum;
  485. extern /* Subroutine */ void slascl_(char *, integer *, integer *, real *,
  486. real *, integer *, integer *, real *, integer *, integer *);
  487. extern integer ilaenv_(integer *, char *, char *, integer *, integer *,
  488. integer *, integer *, ftnlen, ftnlen);
  489. integer ijobvl, iright;
  490. extern /* Subroutine */ void sgeqrf_(integer *, integer *, real *, integer
  491. *, real *, real *, integer *, integer *);
  492. integer ijobvr;
  493. extern /* Subroutine */ void slacpy_(char *, integer *, integer *, real *,
  494. integer *, real *, integer *), slaset_(char *, integer *,
  495. integer *, real *, real *, real *, integer *), stgevc_(
  496. char *, char *, logical *, integer *, real *, integer *, real *,
  497. integer *, real *, integer *, real *, integer *, integer *,
  498. integer *, real *, integer *);
  499. real anrmto, bnrmto;
  500. extern /* Subroutine */ void shgeqz_(char *, char *, char *, integer *,
  501. integer *, integer *, real *, integer *, real *, integer *, real *
  502. , real *, real *, real *, integer *, real *, integer *, real *,
  503. integer *, integer *);
  504. integer minwrk, maxwrk;
  505. real smlnum;
  506. extern /* Subroutine */ void sorgqr_(integer *, integer *, integer *, real
  507. *, integer *, real *, real *, integer *, integer *);
  508. logical lquery;
  509. extern /* Subroutine */ void sormqr_(char *, char *, integer *, integer *,
  510. integer *, real *, integer *, real *, real *, integer *, real *,
  511. integer *, integer *);
  512. integer ihi, ilo;
  513. real eps;
  514. logical ilv;
  515. /* -- LAPACK driver routine (version 3.7.0) -- */
  516. /* -- LAPACK is a software package provided by Univ. of Tennessee, -- */
  517. /* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
  518. /* April 2012 */
  519. /* ===================================================================== */
  520. /* Decode the input arguments */
  521. /* Parameter adjustments */
  522. a_dim1 = *lda;
  523. a_offset = 1 + a_dim1 * 1;
  524. a -= a_offset;
  525. b_dim1 = *ldb;
  526. b_offset = 1 + b_dim1 * 1;
  527. b -= b_offset;
  528. --alphar;
  529. --alphai;
  530. --beta;
  531. vl_dim1 = *ldvl;
  532. vl_offset = 1 + vl_dim1 * 1;
  533. vl -= vl_offset;
  534. vr_dim1 = *ldvr;
  535. vr_offset = 1 + vr_dim1 * 1;
  536. vr -= vr_offset;
  537. --work;
  538. /* Function Body */
  539. if (lsame_(jobvl, "N")) {
  540. ijobvl = 1;
  541. ilvl = FALSE_;
  542. } else if (lsame_(jobvl, "V")) {
  543. ijobvl = 2;
  544. ilvl = TRUE_;
  545. } else {
  546. ijobvl = -1;
  547. ilvl = FALSE_;
  548. }
  549. if (lsame_(jobvr, "N")) {
  550. ijobvr = 1;
  551. ilvr = FALSE_;
  552. } else if (lsame_(jobvr, "V")) {
  553. ijobvr = 2;
  554. ilvr = TRUE_;
  555. } else {
  556. ijobvr = -1;
  557. ilvr = FALSE_;
  558. }
  559. ilv = ilvl || ilvr;
  560. /* Test the input arguments */
  561. *info = 0;
  562. lquery = *lwork == -1;
  563. if (ijobvl <= 0) {
  564. *info = -1;
  565. } else if (ijobvr <= 0) {
  566. *info = -2;
  567. } else if (*n < 0) {
  568. *info = -3;
  569. } else if (*lda < f2cmax(1,*n)) {
  570. *info = -5;
  571. } else if (*ldb < f2cmax(1,*n)) {
  572. *info = -7;
  573. } else if (*ldvl < 1 || ilvl && *ldvl < *n) {
  574. *info = -12;
  575. } else if (*ldvr < 1 || ilvr && *ldvr < *n) {
  576. *info = -14;
  577. }
  578. /* Compute workspace */
  579. /* (Note: Comments in the code beginning "Workspace:" describe the */
  580. /* minimal amount of workspace needed at that point in the code, */
  581. /* as well as the preferred amount for good performance. */
  582. /* NB refers to the optimal block size for the immediately */
  583. /* following subroutine, as returned by ILAENV. The workspace is */
  584. /* computed assuming ILO = 1 and IHI = N, the worst case.) */
  585. if (*info == 0) {
  586. /* Computing MAX */
  587. i__1 = 1, i__2 = *n << 3;
  588. minwrk = f2cmax(i__1,i__2);
  589. /* Computing MAX */
  590. i__1 = 1, i__2 = *n * (ilaenv_(&c__1, "SGEQRF", " ", n, &c__1, n, &
  591. c__0, (ftnlen)6, (ftnlen)1) + 7);
  592. maxwrk = f2cmax(i__1,i__2);
  593. /* Computing MAX */
  594. i__1 = maxwrk, i__2 = *n * (ilaenv_(&c__1, "SORMQR", " ", n, &c__1, n,
  595. &c__0, (ftnlen)6, (ftnlen)1) + 7);
  596. maxwrk = f2cmax(i__1,i__2);
  597. if (ilvl) {
  598. /* Computing MAX */
  599. i__1 = maxwrk, i__2 = *n * (ilaenv_(&c__1, "SORGQR", " ", n, &
  600. c__1, n, &c_n1, (ftnlen)6, (ftnlen)1) + 7);
  601. maxwrk = f2cmax(i__1,i__2);
  602. }
  603. work[1] = (real) maxwrk;
  604. if (*lwork < minwrk && ! lquery) {
  605. *info = -16;
  606. }
  607. }
  608. if (*info != 0) {
  609. i__1 = -(*info);
  610. xerbla_("SGGEV ", &i__1, (ftnlen)6);
  611. return;
  612. } else if (lquery) {
  613. return;
  614. }
  615. /* Quick return if possible */
  616. if (*n == 0) {
  617. return;
  618. }
  619. /* Get machine constants */
  620. eps = slamch_("P");
  621. smlnum = slamch_("S");
  622. bignum = 1.f / smlnum;
  623. slabad_(&smlnum, &bignum);
  624. smlnum = sqrt(smlnum) / eps;
  625. bignum = 1.f / smlnum;
  626. /* Scale A if f2cmax element outside range [SMLNUM,BIGNUM] */
  627. anrm = slange_("M", n, n, &a[a_offset], lda, &work[1]);
  628. ilascl = FALSE_;
  629. if (anrm > 0.f && anrm < smlnum) {
  630. anrmto = smlnum;
  631. ilascl = TRUE_;
  632. } else if (anrm > bignum) {
  633. anrmto = bignum;
  634. ilascl = TRUE_;
  635. }
  636. if (ilascl) {
  637. slascl_("G", &c__0, &c__0, &anrm, &anrmto, n, n, &a[a_offset], lda, &
  638. ierr);
  639. }
  640. /* Scale B if f2cmax element outside range [SMLNUM,BIGNUM] */
  641. bnrm = slange_("M", n, n, &b[b_offset], ldb, &work[1]);
  642. ilbscl = FALSE_;
  643. if (bnrm > 0.f && bnrm < smlnum) {
  644. bnrmto = smlnum;
  645. ilbscl = TRUE_;
  646. } else if (bnrm > bignum) {
  647. bnrmto = bignum;
  648. ilbscl = TRUE_;
  649. }
  650. if (ilbscl) {
  651. slascl_("G", &c__0, &c__0, &bnrm, &bnrmto, n, n, &b[b_offset], ldb, &
  652. ierr);
  653. }
  654. /* Permute the matrices A, B to isolate eigenvalues if possible */
  655. /* (Workspace: need 6*N) */
  656. ileft = 1;
  657. iright = *n + 1;
  658. iwrk = iright + *n;
  659. sggbal_("P", n, &a[a_offset], lda, &b[b_offset], ldb, &ilo, &ihi, &work[
  660. ileft], &work[iright], &work[iwrk], &ierr);
  661. /* Reduce B to triangular form (QR decomposition of B) */
  662. /* (Workspace: need N, prefer N*NB) */
  663. irows = ihi + 1 - ilo;
  664. if (ilv) {
  665. icols = *n + 1 - ilo;
  666. } else {
  667. icols = irows;
  668. }
  669. itau = iwrk;
  670. iwrk = itau + irows;
  671. i__1 = *lwork + 1 - iwrk;
  672. sgeqrf_(&irows, &icols, &b[ilo + ilo * b_dim1], ldb, &work[itau], &work[
  673. iwrk], &i__1, &ierr);
  674. /* Apply the orthogonal transformation to matrix A */
  675. /* (Workspace: need N, prefer N*NB) */
  676. i__1 = *lwork + 1 - iwrk;
  677. sormqr_("L", "T", &irows, &icols, &irows, &b[ilo + ilo * b_dim1], ldb, &
  678. work[itau], &a[ilo + ilo * a_dim1], lda, &work[iwrk], &i__1, &
  679. ierr);
  680. /* Initialize VL */
  681. /* (Workspace: need N, prefer N*NB) */
  682. if (ilvl) {
  683. slaset_("Full", n, n, &c_b36, &c_b37, &vl[vl_offset], ldvl)
  684. ;
  685. if (irows > 1) {
  686. i__1 = irows - 1;
  687. i__2 = irows - 1;
  688. slacpy_("L", &i__1, &i__2, &b[ilo + 1 + ilo * b_dim1], ldb, &vl[
  689. ilo + 1 + ilo * vl_dim1], ldvl);
  690. }
  691. i__1 = *lwork + 1 - iwrk;
  692. sorgqr_(&irows, &irows, &irows, &vl[ilo + ilo * vl_dim1], ldvl, &work[
  693. itau], &work[iwrk], &i__1, &ierr);
  694. }
  695. /* Initialize VR */
  696. if (ilvr) {
  697. slaset_("Full", n, n, &c_b36, &c_b37, &vr[vr_offset], ldvr)
  698. ;
  699. }
  700. /* Reduce to generalized Hessenberg form */
  701. /* (Workspace: none needed) */
  702. if (ilv) {
  703. /* Eigenvectors requested -- work on whole matrix. */
  704. sgghrd_(jobvl, jobvr, n, &ilo, &ihi, &a[a_offset], lda, &b[b_offset],
  705. ldb, &vl[vl_offset], ldvl, &vr[vr_offset], ldvr, &ierr);
  706. } else {
  707. sgghrd_("N", "N", &irows, &c__1, &irows, &a[ilo + ilo * a_dim1], lda,
  708. &b[ilo + ilo * b_dim1], ldb, &vl[vl_offset], ldvl, &vr[
  709. vr_offset], ldvr, &ierr);
  710. }
  711. /* Perform QZ algorithm (Compute eigenvalues, and optionally, the */
  712. /* Schur forms and Schur vectors) */
  713. /* (Workspace: need N) */
  714. iwrk = itau;
  715. if (ilv) {
  716. *(unsigned char *)chtemp = 'S';
  717. } else {
  718. *(unsigned char *)chtemp = 'E';
  719. }
  720. i__1 = *lwork + 1 - iwrk;
  721. shgeqz_(chtemp, jobvl, jobvr, n, &ilo, &ihi, &a[a_offset], lda, &b[
  722. b_offset], ldb, &alphar[1], &alphai[1], &beta[1], &vl[vl_offset],
  723. ldvl, &vr[vr_offset], ldvr, &work[iwrk], &i__1, &ierr);
  724. if (ierr != 0) {
  725. if (ierr > 0 && ierr <= *n) {
  726. *info = ierr;
  727. } else if (ierr > *n && ierr <= *n << 1) {
  728. *info = ierr - *n;
  729. } else {
  730. *info = *n + 1;
  731. }
  732. goto L110;
  733. }
  734. /* Compute Eigenvectors */
  735. /* (Workspace: need 6*N) */
  736. if (ilv) {
  737. if (ilvl) {
  738. if (ilvr) {
  739. *(unsigned char *)chtemp = 'B';
  740. } else {
  741. *(unsigned char *)chtemp = 'L';
  742. }
  743. } else {
  744. *(unsigned char *)chtemp = 'R';
  745. }
  746. stgevc_(chtemp, "B", ldumma, n, &a[a_offset], lda, &b[b_offset], ldb,
  747. &vl[vl_offset], ldvl, &vr[vr_offset], ldvr, n, &in, &work[
  748. iwrk], &ierr);
  749. if (ierr != 0) {
  750. *info = *n + 2;
  751. goto L110;
  752. }
  753. /* Undo balancing on VL and VR and normalization */
  754. /* (Workspace: none needed) */
  755. if (ilvl) {
  756. sggbak_("P", "L", n, &ilo, &ihi, &work[ileft], &work[iright], n, &
  757. vl[vl_offset], ldvl, &ierr);
  758. i__1 = *n;
  759. for (jc = 1; jc <= i__1; ++jc) {
  760. if (alphai[jc] < 0.f) {
  761. goto L50;
  762. }
  763. temp = 0.f;
  764. if (alphai[jc] == 0.f) {
  765. i__2 = *n;
  766. for (jr = 1; jr <= i__2; ++jr) {
  767. /* Computing MAX */
  768. r__2 = temp, r__3 = (r__1 = vl[jr + jc * vl_dim1],
  769. abs(r__1));
  770. temp = f2cmax(r__2,r__3);
  771. /* L10: */
  772. }
  773. } else {
  774. i__2 = *n;
  775. for (jr = 1; jr <= i__2; ++jr) {
  776. /* Computing MAX */
  777. r__3 = temp, r__4 = (r__1 = vl[jr + jc * vl_dim1],
  778. abs(r__1)) + (r__2 = vl[jr + (jc + 1) *
  779. vl_dim1], abs(r__2));
  780. temp = f2cmax(r__3,r__4);
  781. /* L20: */
  782. }
  783. }
  784. if (temp < smlnum) {
  785. goto L50;
  786. }
  787. temp = 1.f / temp;
  788. if (alphai[jc] == 0.f) {
  789. i__2 = *n;
  790. for (jr = 1; jr <= i__2; ++jr) {
  791. vl[jr + jc * vl_dim1] *= temp;
  792. /* L30: */
  793. }
  794. } else {
  795. i__2 = *n;
  796. for (jr = 1; jr <= i__2; ++jr) {
  797. vl[jr + jc * vl_dim1] *= temp;
  798. vl[jr + (jc + 1) * vl_dim1] *= temp;
  799. /* L40: */
  800. }
  801. }
  802. L50:
  803. ;
  804. }
  805. }
  806. if (ilvr) {
  807. sggbak_("P", "R", n, &ilo, &ihi, &work[ileft], &work[iright], n, &
  808. vr[vr_offset], ldvr, &ierr);
  809. i__1 = *n;
  810. for (jc = 1; jc <= i__1; ++jc) {
  811. if (alphai[jc] < 0.f) {
  812. goto L100;
  813. }
  814. temp = 0.f;
  815. if (alphai[jc] == 0.f) {
  816. i__2 = *n;
  817. for (jr = 1; jr <= i__2; ++jr) {
  818. /* Computing MAX */
  819. r__2 = temp, r__3 = (r__1 = vr[jr + jc * vr_dim1],
  820. abs(r__1));
  821. temp = f2cmax(r__2,r__3);
  822. /* L60: */
  823. }
  824. } else {
  825. i__2 = *n;
  826. for (jr = 1; jr <= i__2; ++jr) {
  827. /* Computing MAX */
  828. r__3 = temp, r__4 = (r__1 = vr[jr + jc * vr_dim1],
  829. abs(r__1)) + (r__2 = vr[jr + (jc + 1) *
  830. vr_dim1], abs(r__2));
  831. temp = f2cmax(r__3,r__4);
  832. /* L70: */
  833. }
  834. }
  835. if (temp < smlnum) {
  836. goto L100;
  837. }
  838. temp = 1.f / temp;
  839. if (alphai[jc] == 0.f) {
  840. i__2 = *n;
  841. for (jr = 1; jr <= i__2; ++jr) {
  842. vr[jr + jc * vr_dim1] *= temp;
  843. /* L80: */
  844. }
  845. } else {
  846. i__2 = *n;
  847. for (jr = 1; jr <= i__2; ++jr) {
  848. vr[jr + jc * vr_dim1] *= temp;
  849. vr[jr + (jc + 1) * vr_dim1] *= temp;
  850. /* L90: */
  851. }
  852. }
  853. L100:
  854. ;
  855. }
  856. }
  857. /* End of eigenvector calculation */
  858. }
  859. /* Undo scaling if necessary */
  860. L110:
  861. if (ilascl) {
  862. slascl_("G", &c__0, &c__0, &anrmto, &anrm, n, &c__1, &alphar[1], n, &
  863. ierr);
  864. slascl_("G", &c__0, &c__0, &anrmto, &anrm, n, &c__1, &alphai[1], n, &
  865. ierr);
  866. }
  867. if (ilbscl) {
  868. slascl_("G", &c__0, &c__0, &bnrmto, &bnrm, n, &c__1, &beta[1], n, &
  869. ierr);
  870. }
  871. work[1] = (real) maxwrk;
  872. return;
  873. /* End of SGGEV */
  874. } /* sggev_ */