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.

zgegs.c 30 kB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024
  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. #if defined(_WIN64)
  18. typedef long long BLASLONG;
  19. typedef unsigned long long BLASULONG;
  20. #else
  21. typedef long BLASLONG;
  22. typedef unsigned long BLASULONG;
  23. #endif
  24. #ifdef LAPACK_ILP64
  25. typedef BLASLONG blasint;
  26. #if defined(_WIN64)
  27. #define blasabs(x) llabs(x)
  28. #else
  29. #define blasabs(x) labs(x)
  30. #endif
  31. #else
  32. typedef int blasint;
  33. #define blasabs(x) abs(x)
  34. #endif
  35. typedef blasint integer;
  36. typedef unsigned int uinteger;
  37. typedef char *address;
  38. typedef short int shortint;
  39. typedef float real;
  40. typedef double doublereal;
  41. typedef struct { real r, i; } complex;
  42. typedef struct { doublereal r, i; } doublecomplex;
  43. static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
  44. static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
  45. static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
  46. static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
  47. #define pCf(z) (*_pCf(z))
  48. #define pCd(z) (*_pCd(z))
  49. typedef int logical;
  50. typedef short int shortlogical;
  51. typedef char logical1;
  52. typedef char integer1;
  53. #define TRUE_ (1)
  54. #define FALSE_ (0)
  55. /* Extern is for use with -E */
  56. #ifndef Extern
  57. #define Extern extern
  58. #endif
  59. /* I/O stuff */
  60. typedef int flag;
  61. typedef int ftnlen;
  62. typedef int ftnint;
  63. /*external read, write*/
  64. typedef struct
  65. { flag cierr;
  66. ftnint ciunit;
  67. flag ciend;
  68. char *cifmt;
  69. ftnint cirec;
  70. } cilist;
  71. /*internal read, write*/
  72. typedef struct
  73. { flag icierr;
  74. char *iciunit;
  75. flag iciend;
  76. char *icifmt;
  77. ftnint icirlen;
  78. ftnint icirnum;
  79. } icilist;
  80. /*open*/
  81. typedef struct
  82. { flag oerr;
  83. ftnint ounit;
  84. char *ofnm;
  85. ftnlen ofnmlen;
  86. char *osta;
  87. char *oacc;
  88. char *ofm;
  89. ftnint orl;
  90. char *oblnk;
  91. } olist;
  92. /*close*/
  93. typedef struct
  94. { flag cerr;
  95. ftnint cunit;
  96. char *csta;
  97. } cllist;
  98. /*rewind, backspace, endfile*/
  99. typedef struct
  100. { flag aerr;
  101. ftnint aunit;
  102. } alist;
  103. /* inquire */
  104. typedef struct
  105. { flag inerr;
  106. ftnint inunit;
  107. char *infile;
  108. ftnlen infilen;
  109. ftnint *inex; /*parameters in standard's order*/
  110. ftnint *inopen;
  111. ftnint *innum;
  112. ftnint *innamed;
  113. char *inname;
  114. ftnlen innamlen;
  115. char *inacc;
  116. ftnlen inacclen;
  117. char *inseq;
  118. ftnlen inseqlen;
  119. char *indir;
  120. ftnlen indirlen;
  121. char *infmt;
  122. ftnlen infmtlen;
  123. char *inform;
  124. ftnint informlen;
  125. char *inunf;
  126. ftnlen inunflen;
  127. ftnint *inrecl;
  128. ftnint *innrec;
  129. char *inblank;
  130. ftnlen inblanklen;
  131. } inlist;
  132. #define VOID void
  133. union Multitype { /* for multiple entry points */
  134. integer1 g;
  135. shortint h;
  136. integer i;
  137. /* longint j; */
  138. real r;
  139. doublereal d;
  140. complex c;
  141. doublecomplex z;
  142. };
  143. typedef union Multitype Multitype;
  144. struct Vardesc { /* for Namelist */
  145. char *name;
  146. char *addr;
  147. ftnlen *dims;
  148. int type;
  149. };
  150. typedef struct Vardesc Vardesc;
  151. struct Namelist {
  152. char *name;
  153. Vardesc **vars;
  154. int nvars;
  155. };
  156. typedef struct Namelist Namelist;
  157. #define abs(x) ((x) >= 0 ? (x) : -(x))
  158. #define dabs(x) (fabs(x))
  159. #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
  160. #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
  161. #define dmin(a,b) (f2cmin(a,b))
  162. #define dmax(a,b) (f2cmax(a,b))
  163. #define bit_test(a,b) ((a) >> (b) & 1)
  164. #define bit_clear(a,b) ((a) & ~((uinteger)1 << (b)))
  165. #define bit_set(a,b) ((a) | ((uinteger)1 << (b)))
  166. #define abort_() { sig_die("Fortran abort routine called", 1); }
  167. #define c_abs(z) (cabsf(Cf(z)))
  168. #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
  169. #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
  170. #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
  171. #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
  172. #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
  173. #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
  174. //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
  175. #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
  176. #define d_abs(x) (fabs(*(x)))
  177. #define d_acos(x) (acos(*(x)))
  178. #define d_asin(x) (asin(*(x)))
  179. #define d_atan(x) (atan(*(x)))
  180. #define d_atn2(x, y) (atan2(*(x),*(y)))
  181. #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
  182. #define r_cnjg(R, Z) { pCf(R) = conj(Cf(Z)); }
  183. #define d_cos(x) (cos(*(x)))
  184. #define d_cosh(x) (cosh(*(x)))
  185. #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
  186. #define d_exp(x) (exp(*(x)))
  187. #define d_imag(z) (cimag(Cd(z)))
  188. #define r_imag(z) (cimag(Cf(z)))
  189. #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
  190. #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
  191. #define d_log(x) (log(*(x)))
  192. #define d_mod(x, y) (fmod(*(x), *(y)))
  193. #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
  194. #define d_nint(x) u_nint(*(x))
  195. #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
  196. #define d_sign(a,b) u_sign(*(a),*(b))
  197. #define d_sin(x) (sin(*(x)))
  198. #define d_sinh(x) (sinh(*(x)))
  199. #define d_sqrt(x) (sqrt(*(x)))
  200. #define d_tan(x) (tan(*(x)))
  201. #define d_tanh(x) (tanh(*(x)))
  202. #define i_abs(x) abs(*(x))
  203. #define i_dnnt(x) ((integer)u_nint(*(x)))
  204. #define i_len(s, n) (n)
  205. #define i_nint(x) ((integer)u_nint(*(x)))
  206. #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
  207. #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
  208. #define pow_si(B,E) spow_ui(*(B),*(E))
  209. #define pow_di(B,E) dpow_ui(*(B),*(E))
  210. #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
  211. #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
  212. #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++ = ' '; }
  213. #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
  214. #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]; }
  215. #define sig_die(s, kill) { exit(1); }
  216. #define s_stop(s, n) {exit(0);}
  217. static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
  218. #define z_abs(z) (cabs(Cd(z)))
  219. #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
  220. #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
  221. #define myexit_() break;
  222. #define mycycle() continue;
  223. #define myceiling(w) {ceil(w)}
  224. #define myhuge(w) {HUGE_VAL}
  225. #define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
  226. /* procedure parameter types for -A and -C++ */
  227. #define F2C_proc_par_types 1
  228. #ifdef __cplusplus
  229. typedef logical (*L_fp)(...);
  230. #else
  231. typedef logical (*L_fp)();
  232. #endif
  233. static float spow_ui(float x, integer n) {
  234. float 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 double dpow_ui(double x, integer n) {
  246. double 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 float cpow_ui(_Complex float x, integer n) {
  258. _Complex float 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 _Complex double zpow_ui(_Complex double x, integer n) {
  270. _Complex double pow=1.0; unsigned long int u;
  271. if(n != 0) {
  272. if(n < 0) n = -n, x = 1/x;
  273. for(u = n; ; ) {
  274. if(u & 01) pow *= x;
  275. if(u >>= 1) x *= x;
  276. else break;
  277. }
  278. }
  279. return pow;
  280. }
  281. static integer pow_ii(integer x, integer n) {
  282. integer pow; unsigned long int u;
  283. if (n <= 0) {
  284. if (n == 0 || x == 1) pow = 1;
  285. else if (x != -1) pow = x == 0 ? 1/x : 0;
  286. else n = -n;
  287. }
  288. if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
  289. u = n;
  290. for(pow = 1; ; ) {
  291. if(u & 01) pow *= x;
  292. if(u >>= 1) x *= x;
  293. else break;
  294. }
  295. }
  296. return pow;
  297. }
  298. static integer dmaxloc_(double *w, integer s, integer e, integer *n)
  299. {
  300. double m; integer i, mi;
  301. for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
  302. if (w[i-1]>m) mi=i ,m=w[i-1];
  303. return mi-s+1;
  304. }
  305. static integer smaxloc_(float *w, integer s, integer e, integer *n)
  306. {
  307. float m; integer i, mi;
  308. for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
  309. if (w[i-1]>m) mi=i ,m=w[i-1];
  310. return mi-s+1;
  311. }
  312. static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
  313. integer n = *n_, incx = *incx_, incy = *incy_, i;
  314. _Complex float zdotc = 0.0;
  315. if (incx == 1 && incy == 1) {
  316. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  317. zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
  318. }
  319. } else {
  320. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  321. zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
  322. }
  323. }
  324. pCf(z) = zdotc;
  325. }
  326. static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
  327. integer n = *n_, incx = *incx_, incy = *incy_, i;
  328. _Complex double zdotc = 0.0;
  329. if (incx == 1 && incy == 1) {
  330. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  331. zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
  332. }
  333. } else {
  334. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  335. zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
  336. }
  337. }
  338. pCd(z) = zdotc;
  339. }
  340. static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
  341. integer n = *n_, incx = *incx_, incy = *incy_, i;
  342. _Complex float zdotc = 0.0;
  343. if (incx == 1 && incy == 1) {
  344. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  345. zdotc += Cf(&x[i]) * Cf(&y[i]);
  346. }
  347. } else {
  348. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  349. zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
  350. }
  351. }
  352. pCf(z) = zdotc;
  353. }
  354. static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
  355. integer n = *n_, incx = *incx_, incy = *incy_, i;
  356. _Complex double zdotc = 0.0;
  357. if (incx == 1 && incy == 1) {
  358. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  359. zdotc += Cd(&x[i]) * Cd(&y[i]);
  360. }
  361. } else {
  362. for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
  363. zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
  364. }
  365. }
  366. pCd(z) = zdotc;
  367. }
  368. #endif
  369. /* -- translated by f2c (version 20000121).
  370. You must link the resulting object file with the libraries:
  371. -lf2c -lm (in that order)
  372. */
  373. /* Table of constant values */
  374. static doublecomplex c_b1 = {0.,0.};
  375. static doublecomplex c_b2 = {1.,0.};
  376. static integer c__1 = 1;
  377. static integer c_n1 = -1;
  378. /* > \brief <b> ZGEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE mat
  379. rices</b> */
  380. /* =========== DOCUMENTATION =========== */
  381. /* Online html documentation available at */
  382. /* http://www.netlib.org/lapack/explore-html/ */
  383. /* > \htmlonly */
  384. /* > Download ZGEGS + dependencies */
  385. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zgegs.f
  386. "> */
  387. /* > [TGZ]</a> */
  388. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zgegs.f
  389. "> */
  390. /* > [ZIP]</a> */
  391. /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zgegs.f
  392. "> */
  393. /* > [TXT]</a> */
  394. /* > \endhtmlonly */
  395. /* Definition: */
  396. /* =========== */
  397. /* SUBROUTINE ZGEGS( JOBVSL, JOBVSR, N, A, LDA, B, LDB, ALPHA, BETA, */
  398. /* VSL, LDVSL, VSR, LDVSR, WORK, LWORK, RWORK, */
  399. /* INFO ) */
  400. /* CHARACTER JOBVSL, JOBVSR */
  401. /* INTEGER INFO, LDA, LDB, LDVSL, LDVSR, LWORK, N */
  402. /* DOUBLE PRECISION RWORK( * ) */
  403. /* COMPLEX*16 A( LDA, * ), ALPHA( * ), B( LDB, * ), */
  404. /* $ BETA( * ), VSL( LDVSL, * ), VSR( LDVSR, * ), */
  405. /* $ WORK( * ) */
  406. /* > \par Purpose: */
  407. /* ============= */
  408. /* > */
  409. /* > \verbatim */
  410. /* > */
  411. /* > This routine is deprecated and has been replaced by routine ZGGES. */
  412. /* > */
  413. /* > ZGEGS computes the eigenvalues, Schur form, and, optionally, the */
  414. /* > left and or/right Schur vectors of a complex matrix pair (A,B). */
  415. /* > Given two square matrices A and B, the generalized Schur */
  416. /* > factorization has the form */
  417. /* > */
  418. /* > A = Q*S*Z**H, B = Q*T*Z**H */
  419. /* > */
  420. /* > where Q and Z are unitary matrices and S and T are upper triangular. */
  421. /* > The columns of Q are the left Schur vectors */
  422. /* > and the columns of Z are the right Schur vectors. */
  423. /* > */
  424. /* > If only the eigenvalues of (A,B) are needed, the driver routine */
  425. /* > ZGEGV should be used instead. See ZGEGV for a description of the */
  426. /* > eigenvalues of the generalized nonsymmetric eigenvalue problem */
  427. /* > (GNEP). */
  428. /* > \endverbatim */
  429. /* Arguments: */
  430. /* ========== */
  431. /* > \param[in] JOBVSL */
  432. /* > \verbatim */
  433. /* > JOBVSL is CHARACTER*1 */
  434. /* > = 'N': do not compute the left Schur vectors; */
  435. /* > = 'V': compute the left Schur vectors (returned in VSL). */
  436. /* > \endverbatim */
  437. /* > */
  438. /* > \param[in] JOBVSR */
  439. /* > \verbatim */
  440. /* > JOBVSR is CHARACTER*1 */
  441. /* > = 'N': do not compute the right Schur vectors; */
  442. /* > = 'V': compute the right Schur vectors (returned in VSR). */
  443. /* > \endverbatim */
  444. /* > */
  445. /* > \param[in] N */
  446. /* > \verbatim */
  447. /* > N is INTEGER */
  448. /* > The order of the matrices A, B, VSL, and VSR. N >= 0. */
  449. /* > \endverbatim */
  450. /* > */
  451. /* > \param[in,out] A */
  452. /* > \verbatim */
  453. /* > A is COMPLEX*16 array, dimension (LDA, N) */
  454. /* > On entry, the matrix A. */
  455. /* > On exit, the upper triangular matrix S from the generalized */
  456. /* > Schur factorization. */
  457. /* > \endverbatim */
  458. /* > */
  459. /* > \param[in] LDA */
  460. /* > \verbatim */
  461. /* > LDA is INTEGER */
  462. /* > The leading dimension of A. LDA >= f2cmax(1,N). */
  463. /* > \endverbatim */
  464. /* > */
  465. /* > \param[in,out] B */
  466. /* > \verbatim */
  467. /* > B is COMPLEX*16 array, dimension (LDB, N) */
  468. /* > On entry, the matrix B. */
  469. /* > On exit, the upper triangular matrix T from the generalized */
  470. /* > Schur factorization. */
  471. /* > \endverbatim */
  472. /* > */
  473. /* > \param[in] LDB */
  474. /* > \verbatim */
  475. /* > LDB is INTEGER */
  476. /* > The leading dimension of B. LDB >= f2cmax(1,N). */
  477. /* > \endverbatim */
  478. /* > */
  479. /* > \param[out] ALPHA */
  480. /* > \verbatim */
  481. /* > ALPHA is COMPLEX*16 array, dimension (N) */
  482. /* > The complex scalars alpha that define the eigenvalues of */
  483. /* > GNEP. ALPHA(j) = S(j,j), the diagonal element of the Schur */
  484. /* > form of A. */
  485. /* > \endverbatim */
  486. /* > */
  487. /* > \param[out] BETA */
  488. /* > \verbatim */
  489. /* > BETA is COMPLEX*16 array, dimension (N) */
  490. /* > The non-negative real scalars beta that define the */
  491. /* > eigenvalues of GNEP. BETA(j) = T(j,j), the diagonal element */
  492. /* > of the triangular factor T. */
  493. /* > */
  494. /* > Together, the quantities alpha = ALPHA(j) and beta = BETA(j) */
  495. /* > represent the j-th eigenvalue of the matrix pair (A,B), in */
  496. /* > one of the forms lambda = alpha/beta or mu = beta/alpha. */
  497. /* > Since either lambda or mu may overflow, they should not, */
  498. /* > in general, be computed. */
  499. /* > \endverbatim */
  500. /* > */
  501. /* > \param[out] VSL */
  502. /* > \verbatim */
  503. /* > VSL is COMPLEX*16 array, dimension (LDVSL,N) */
  504. /* > If JOBVSL = 'V', the matrix of left Schur vectors Q. */
  505. /* > Not referenced if JOBVSL = 'N'. */
  506. /* > \endverbatim */
  507. /* > */
  508. /* > \param[in] LDVSL */
  509. /* > \verbatim */
  510. /* > LDVSL is INTEGER */
  511. /* > The leading dimension of the matrix VSL. LDVSL >= 1, and */
  512. /* > if JOBVSL = 'V', LDVSL >= N. */
  513. /* > \endverbatim */
  514. /* > */
  515. /* > \param[out] VSR */
  516. /* > \verbatim */
  517. /* > VSR is COMPLEX*16 array, dimension (LDVSR,N) */
  518. /* > If JOBVSR = 'V', the matrix of right Schur vectors Z. */
  519. /* > Not referenced if JOBVSR = 'N'. */
  520. /* > \endverbatim */
  521. /* > */
  522. /* > \param[in] LDVSR */
  523. /* > \verbatim */
  524. /* > LDVSR is INTEGER */
  525. /* > The leading dimension of the matrix VSR. LDVSR >= 1, and */
  526. /* > if JOBVSR = 'V', LDVSR >= N. */
  527. /* > \endverbatim */
  528. /* > */
  529. /* > \param[out] WORK */
  530. /* > \verbatim */
  531. /* > WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)) */
  532. /* > On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
  533. /* > \endverbatim */
  534. /* > */
  535. /* > \param[in] LWORK */
  536. /* > \verbatim */
  537. /* > LWORK is INTEGER */
  538. /* > The dimension of the array WORK. LWORK >= f2cmax(1,2*N). */
  539. /* > For good performance, LWORK must generally be larger. */
  540. /* > To compute the optimal value of LWORK, call ILAENV to get */
  541. /* > blocksizes (for ZGEQRF, ZUNMQR, and CUNGQR.) Then compute: */
  542. /* > NB -- MAX of the blocksizes for ZGEQRF, ZUNMQR, and CUNGQR; */
  543. /* > the optimal LWORK is N*(NB+1). */
  544. /* > */
  545. /* > If LWORK = -1, then a workspace query is assumed; the routine */
  546. /* > only calculates the optimal size of the WORK array, returns */
  547. /* > this value as the first entry of the WORK array, and no error */
  548. /* > message related to LWORK is issued by XERBLA. */
  549. /* > \endverbatim */
  550. /* > */
  551. /* > \param[out] RWORK */
  552. /* > \verbatim */
  553. /* > RWORK is DOUBLE PRECISION array, dimension (3*N) */
  554. /* > \endverbatim */
  555. /* > */
  556. /* > \param[out] INFO */
  557. /* > \verbatim */
  558. /* > INFO is INTEGER */
  559. /* > = 0: successful exit */
  560. /* > < 0: if INFO = -i, the i-th argument had an illegal value. */
  561. /* > =1,...,N: */
  562. /* > The QZ iteration failed. (A,B) are not in Schur */
  563. /* > form, but ALPHA(j) and BETA(j) should be correct for */
  564. /* > j=INFO+1,...,N. */
  565. /* > > N: errors that usually indicate LAPACK problems: */
  566. /* > =N+1: error return from ZGGBAL */
  567. /* > =N+2: error return from ZGEQRF */
  568. /* > =N+3: error return from ZUNMQR */
  569. /* > =N+4: error return from ZUNGQR */
  570. /* > =N+5: error return from ZGGHRD */
  571. /* > =N+6: error return from ZHGEQZ (other than failed */
  572. /* > iteration) */
  573. /* > =N+7: error return from ZGGBAK (computing VSL) */
  574. /* > =N+8: error return from ZGGBAK (computing VSR) */
  575. /* > =N+9: error return from ZLASCL (various places) */
  576. /* > \endverbatim */
  577. /* Authors: */
  578. /* ======== */
  579. /* > \author Univ. of Tennessee */
  580. /* > \author Univ. of California Berkeley */
  581. /* > \author Univ. of Colorado Denver */
  582. /* > \author NAG Ltd. */
  583. /* > \date December 2016 */
  584. /* > \ingroup complex16GEeigen */
  585. /* ===================================================================== */
  586. /* Subroutine */ int zgegs_(char *jobvsl, char *jobvsr, integer *n,
  587. doublecomplex *a, integer *lda, doublecomplex *b, integer *ldb,
  588. doublecomplex *alpha, doublecomplex *beta, doublecomplex *vsl,
  589. integer *ldvsl, doublecomplex *vsr, integer *ldvsr, doublecomplex *
  590. work, integer *lwork, doublereal *rwork, integer *info)
  591. {
  592. /* System generated locals */
  593. integer a_dim1, a_offset, b_dim1, b_offset, vsl_dim1, vsl_offset,
  594. vsr_dim1, vsr_offset, i__1, i__2, i__3;
  595. /* Local variables */
  596. doublereal anrm, bnrm;
  597. integer itau, lopt;
  598. extern logical lsame_(char *, char *);
  599. integer ileft, iinfo, icols;
  600. logical ilvsl;
  601. integer iwork;
  602. logical ilvsr;
  603. integer irows, nb;
  604. extern doublereal dlamch_(char *);
  605. extern /* Subroutine */ int zggbak_(char *, char *, integer *, integer *,
  606. integer *, doublereal *, doublereal *, integer *, doublecomplex *,
  607. integer *, integer *), zggbal_(char *, integer *,
  608. doublecomplex *, integer *, doublecomplex *, integer *, integer *
  609. , integer *, doublereal *, doublereal *, doublereal *, integer *);
  610. logical ilascl, ilbscl;
  611. doublereal safmin;
  612. extern /* Subroutine */ int xerbla_(char *, integer *);
  613. extern integer ilaenv_(integer *, char *, char *, integer *, integer *,
  614. integer *, integer *, ftnlen, ftnlen);
  615. extern doublereal zlange_(char *, integer *, integer *, doublecomplex *,
  616. integer *, doublereal *);
  617. doublereal bignum;
  618. integer ijobvl, iright;
  619. extern /* Subroutine */ int zgghrd_(char *, char *, integer *, integer *,
  620. integer *, doublecomplex *, integer *, doublecomplex *, integer *,
  621. doublecomplex *, integer *, doublecomplex *, integer *, integer *
  622. ), zlascl_(char *, integer *, integer *,
  623. doublereal *, doublereal *, integer *, integer *, doublecomplex *,
  624. integer *, integer *);
  625. integer ijobvr;
  626. extern /* Subroutine */ int zgeqrf_(integer *, integer *, doublecomplex *,
  627. integer *, doublecomplex *, doublecomplex *, integer *, integer *
  628. );
  629. doublereal anrmto;
  630. integer lwkmin, nb1, nb2, nb3;
  631. doublereal bnrmto;
  632. extern /* Subroutine */ int zlacpy_(char *, integer *, integer *,
  633. doublecomplex *, integer *, doublecomplex *, integer *),
  634. zhgeqz_(char *, char *, char *, integer *, integer *, integer *,
  635. doublecomplex *, integer *, doublecomplex *, integer *,
  636. doublecomplex *, doublecomplex *, doublecomplex *, integer *,
  637. doublecomplex *, integer *, doublecomplex *, integer *,
  638. doublereal *, integer *), zlaset_(char *,
  639. integer *, integer *, doublecomplex *, doublecomplex *,
  640. doublecomplex *, integer *);
  641. doublereal smlnum;
  642. integer irwork, lwkopt;
  643. logical lquery;
  644. extern /* Subroutine */ int zungqr_(integer *, integer *, integer *,
  645. doublecomplex *, integer *, doublecomplex *, doublecomplex *,
  646. integer *, integer *), zunmqr_(char *, char *, integer *, integer
  647. *, integer *, doublecomplex *, integer *, doublecomplex *,
  648. doublecomplex *, integer *, doublecomplex *, integer *, integer *);
  649. integer ihi, ilo;
  650. doublereal eps;
  651. /* -- LAPACK driver routine (version 3.7.0) -- */
  652. /* -- LAPACK is a software package provided by Univ. of Tennessee, -- */
  653. /* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
  654. /* December 2016 */
  655. /* ===================================================================== */
  656. /* Decode the input arguments */
  657. /* Parameter adjustments */
  658. a_dim1 = *lda;
  659. a_offset = 1 + a_dim1 * 1;
  660. a -= a_offset;
  661. b_dim1 = *ldb;
  662. b_offset = 1 + b_dim1 * 1;
  663. b -= b_offset;
  664. --alpha;
  665. --beta;
  666. vsl_dim1 = *ldvsl;
  667. vsl_offset = 1 + vsl_dim1 * 1;
  668. vsl -= vsl_offset;
  669. vsr_dim1 = *ldvsr;
  670. vsr_offset = 1 + vsr_dim1 * 1;
  671. vsr -= vsr_offset;
  672. --work;
  673. --rwork;
  674. /* Function Body */
  675. if (lsame_(jobvsl, "N")) {
  676. ijobvl = 1;
  677. ilvsl = FALSE_;
  678. } else if (lsame_(jobvsl, "V")) {
  679. ijobvl = 2;
  680. ilvsl = TRUE_;
  681. } else {
  682. ijobvl = -1;
  683. ilvsl = FALSE_;
  684. }
  685. if (lsame_(jobvsr, "N")) {
  686. ijobvr = 1;
  687. ilvsr = FALSE_;
  688. } else if (lsame_(jobvsr, "V")) {
  689. ijobvr = 2;
  690. ilvsr = TRUE_;
  691. } else {
  692. ijobvr = -1;
  693. ilvsr = FALSE_;
  694. }
  695. /* Test the input arguments */
  696. /* Computing MAX */
  697. i__1 = *n << 1;
  698. lwkmin = f2cmax(i__1,1);
  699. lwkopt = lwkmin;
  700. work[1].r = (doublereal) lwkopt, work[1].i = 0.;
  701. lquery = *lwork == -1;
  702. *info = 0;
  703. if (ijobvl <= 0) {
  704. *info = -1;
  705. } else if (ijobvr <= 0) {
  706. *info = -2;
  707. } else if (*n < 0) {
  708. *info = -3;
  709. } else if (*lda < f2cmax(1,*n)) {
  710. *info = -5;
  711. } else if (*ldb < f2cmax(1,*n)) {
  712. *info = -7;
  713. } else if (*ldvsl < 1 || ilvsl && *ldvsl < *n) {
  714. *info = -11;
  715. } else if (*ldvsr < 1 || ilvsr && *ldvsr < *n) {
  716. *info = -13;
  717. } else if (*lwork < lwkmin && ! lquery) {
  718. *info = -15;
  719. }
  720. if (*info == 0) {
  721. nb1 = ilaenv_(&c__1, "ZGEQRF", " ", n, n, &c_n1, &c_n1, (ftnlen)6, (
  722. ftnlen)1);
  723. nb2 = ilaenv_(&c__1, "ZUNMQR", " ", n, n, n, &c_n1, (ftnlen)6, (
  724. ftnlen)1);
  725. nb3 = ilaenv_(&c__1, "ZUNGQR", " ", n, n, n, &c_n1, (ftnlen)6, (
  726. ftnlen)1);
  727. /* Computing MAX */
  728. i__1 = f2cmax(nb1,nb2);
  729. nb = f2cmax(i__1,nb3);
  730. lopt = *n * (nb + 1);
  731. work[1].r = (doublereal) lopt, work[1].i = 0.;
  732. }
  733. if (*info != 0) {
  734. i__1 = -(*info);
  735. xerbla_("ZGEGS ", &i__1);
  736. return 0;
  737. } else if (lquery) {
  738. return 0;
  739. }
  740. /* Quick return if possible */
  741. if (*n == 0) {
  742. return 0;
  743. }
  744. /* Get machine constants */
  745. eps = dlamch_("E") * dlamch_("B");
  746. safmin = dlamch_("S");
  747. smlnum = *n * safmin / eps;
  748. bignum = 1. / smlnum;
  749. /* Scale A if f2cmax element outside range [SMLNUM,BIGNUM] */
  750. anrm = zlange_("M", n, n, &a[a_offset], lda, &rwork[1]);
  751. ilascl = FALSE_;
  752. if (anrm > 0. && anrm < smlnum) {
  753. anrmto = smlnum;
  754. ilascl = TRUE_;
  755. } else if (anrm > bignum) {
  756. anrmto = bignum;
  757. ilascl = TRUE_;
  758. }
  759. if (ilascl) {
  760. zlascl_("G", &c_n1, &c_n1, &anrm, &anrmto, n, n, &a[a_offset], lda, &
  761. iinfo);
  762. if (iinfo != 0) {
  763. *info = *n + 9;
  764. return 0;
  765. }
  766. }
  767. /* Scale B if f2cmax element outside range [SMLNUM,BIGNUM] */
  768. bnrm = zlange_("M", n, n, &b[b_offset], ldb, &rwork[1]);
  769. ilbscl = FALSE_;
  770. if (bnrm > 0. && bnrm < smlnum) {
  771. bnrmto = smlnum;
  772. ilbscl = TRUE_;
  773. } else if (bnrm > bignum) {
  774. bnrmto = bignum;
  775. ilbscl = TRUE_;
  776. }
  777. if (ilbscl) {
  778. zlascl_("G", &c_n1, &c_n1, &bnrm, &bnrmto, n, n, &b[b_offset], ldb, &
  779. iinfo);
  780. if (iinfo != 0) {
  781. *info = *n + 9;
  782. return 0;
  783. }
  784. }
  785. /* Permute the matrix to make it more nearly triangular */
  786. ileft = 1;
  787. iright = *n + 1;
  788. irwork = iright + *n;
  789. iwork = 1;
  790. zggbal_("P", n, &a[a_offset], lda, &b[b_offset], ldb, &ilo, &ihi, &rwork[
  791. ileft], &rwork[iright], &rwork[irwork], &iinfo);
  792. if (iinfo != 0) {
  793. *info = *n + 1;
  794. goto L10;
  795. }
  796. /* Reduce B to triangular form, and initialize VSL and/or VSR */
  797. irows = ihi + 1 - ilo;
  798. icols = *n + 1 - ilo;
  799. itau = iwork;
  800. iwork = itau + irows;
  801. i__1 = *lwork + 1 - iwork;
  802. zgeqrf_(&irows, &icols, &b[ilo + ilo * b_dim1], ldb, &work[itau], &work[
  803. iwork], &i__1, &iinfo);
  804. if (iinfo >= 0) {
  805. /* Computing MAX */
  806. i__3 = iwork;
  807. i__1 = lwkopt, i__2 = (integer) work[i__3].r + iwork - 1;
  808. lwkopt = f2cmax(i__1,i__2);
  809. }
  810. if (iinfo != 0) {
  811. *info = *n + 2;
  812. goto L10;
  813. }
  814. i__1 = *lwork + 1 - iwork;
  815. zunmqr_("L", "C", &irows, &icols, &irows, &b[ilo + ilo * b_dim1], ldb, &
  816. work[itau], &a[ilo + ilo * a_dim1], lda, &work[iwork], &i__1, &
  817. iinfo);
  818. if (iinfo >= 0) {
  819. /* Computing MAX */
  820. i__3 = iwork;
  821. i__1 = lwkopt, i__2 = (integer) work[i__3].r + iwork - 1;
  822. lwkopt = f2cmax(i__1,i__2);
  823. }
  824. if (iinfo != 0) {
  825. *info = *n + 3;
  826. goto L10;
  827. }
  828. if (ilvsl) {
  829. zlaset_("Full", n, n, &c_b1, &c_b2, &vsl[vsl_offset], ldvsl);
  830. i__1 = irows - 1;
  831. i__2 = irows - 1;
  832. zlacpy_("L", &i__1, &i__2, &b[ilo + 1 + ilo * b_dim1], ldb, &vsl[ilo
  833. + 1 + ilo * vsl_dim1], ldvsl);
  834. i__1 = *lwork + 1 - iwork;
  835. zungqr_(&irows, &irows, &irows, &vsl[ilo + ilo * vsl_dim1], ldvsl, &
  836. work[itau], &work[iwork], &i__1, &iinfo);
  837. if (iinfo >= 0) {
  838. /* Computing MAX */
  839. i__3 = iwork;
  840. i__1 = lwkopt, i__2 = (integer) work[i__3].r + iwork - 1;
  841. lwkopt = f2cmax(i__1,i__2);
  842. }
  843. if (iinfo != 0) {
  844. *info = *n + 4;
  845. goto L10;
  846. }
  847. }
  848. if (ilvsr) {
  849. zlaset_("Full", n, n, &c_b1, &c_b2, &vsr[vsr_offset], ldvsr);
  850. }
  851. /* Reduce to generalized Hessenberg form */
  852. zgghrd_(jobvsl, jobvsr, n, &ilo, &ihi, &a[a_offset], lda, &b[b_offset],
  853. ldb, &vsl[vsl_offset], ldvsl, &vsr[vsr_offset], ldvsr, &iinfo);
  854. if (iinfo != 0) {
  855. *info = *n + 5;
  856. goto L10;
  857. }
  858. /* Perform QZ algorithm, computing Schur vectors if desired */
  859. iwork = itau;
  860. i__1 = *lwork + 1 - iwork;
  861. zhgeqz_("S", jobvsl, jobvsr, n, &ilo, &ihi, &a[a_offset], lda, &b[
  862. b_offset], ldb, &alpha[1], &beta[1], &vsl[vsl_offset], ldvsl, &
  863. vsr[vsr_offset], ldvsr, &work[iwork], &i__1, &rwork[irwork], &
  864. iinfo);
  865. if (iinfo >= 0) {
  866. /* Computing MAX */
  867. i__3 = iwork;
  868. i__1 = lwkopt, i__2 = (integer) work[i__3].r + iwork - 1;
  869. lwkopt = f2cmax(i__1,i__2);
  870. }
  871. if (iinfo != 0) {
  872. if (iinfo > 0 && iinfo <= *n) {
  873. *info = iinfo;
  874. } else if (iinfo > *n && iinfo <= *n << 1) {
  875. *info = iinfo - *n;
  876. } else {
  877. *info = *n + 6;
  878. }
  879. goto L10;
  880. }
  881. /* Apply permutation to VSL and VSR */
  882. if (ilvsl) {
  883. zggbak_("P", "L", n, &ilo, &ihi, &rwork[ileft], &rwork[iright], n, &
  884. vsl[vsl_offset], ldvsl, &iinfo);
  885. if (iinfo != 0) {
  886. *info = *n + 7;
  887. goto L10;
  888. }
  889. }
  890. if (ilvsr) {
  891. zggbak_("P", "R", n, &ilo, &ihi, &rwork[ileft], &rwork[iright], n, &
  892. vsr[vsr_offset], ldvsr, &iinfo);
  893. if (iinfo != 0) {
  894. *info = *n + 8;
  895. goto L10;
  896. }
  897. }
  898. /* Undo scaling */
  899. if (ilascl) {
  900. zlascl_("U", &c_n1, &c_n1, &anrmto, &anrm, n, n, &a[a_offset], lda, &
  901. iinfo);
  902. if (iinfo != 0) {
  903. *info = *n + 9;
  904. return 0;
  905. }
  906. zlascl_("G", &c_n1, &c_n1, &anrmto, &anrm, n, &c__1, &alpha[1], n, &
  907. iinfo);
  908. if (iinfo != 0) {
  909. *info = *n + 9;
  910. return 0;
  911. }
  912. }
  913. if (ilbscl) {
  914. zlascl_("U", &c_n1, &c_n1, &bnrmto, &bnrm, n, n, &b[b_offset], ldb, &
  915. iinfo);
  916. if (iinfo != 0) {
  917. *info = *n + 9;
  918. return 0;
  919. }
  920. zlascl_("G", &c_n1, &c_n1, &bnrmto, &bnrm, n, &c__1, &beta[1], n, &
  921. iinfo);
  922. if (iinfo != 0) {
  923. *info = *n + 9;
  924. return 0;
  925. }
  926. }
  927. L10:
  928. work[1].r = (doublereal) lwkopt, work[1].i = 0.;
  929. return 0;
  930. /* End of ZGEGS */
  931. } /* zgegs_ */