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 26 kB

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