14 typedef long long BLASLONG;
15 typedef unsigned long long BLASULONG;
17 typedef long BLASLONG;
18 typedef unsigned long BLASULONG;
22 typedef BLASLONG blasint;
24 #define blasabs(x) llabs(x)
26 #define blasabs(x) labs(x)
30 #define blasabs(x) abs(x)
33 typedef blasint integer;
35 typedef unsigned int uinteger;
36 typedef char *address;
37 typedef short int shortint;
39 typedef double doublereal;
40 typedef struct { real r, i; } complex;
41 typedef struct { doublereal r, i; } doublecomplex;
43 static inline _Fcomplex Cf(complex *z) {_Fcomplex zz={z->r , z->i}; return zz;}
44 static inline _Dcomplex Cd(doublecomplex *z) {_Dcomplex zz={z->r , z->i};return zz;}
45 static inline _Fcomplex * _pCf(complex *z) {return (_Fcomplex*)z;}
46 static inline _Dcomplex * _pCd(doublecomplex *z) {return (_Dcomplex*)z;}
48 static inline _Complex float Cf(complex *z) {return z->r + z->i*_Complex_I;}
49 static inline _Complex double Cd(doublecomplex *z) {return z->r + z->i*_Complex_I;}
50 static inline _Complex float * _pCf(complex *z) {return (_Complex float*)z;}
51 static inline _Complex double * _pCd(doublecomplex *z) {return (_Complex double*)z;}
53 #define pCf(z) (*_pCf(z))
54 #define pCd(z) (*_pCd(z))
56 typedef short int shortlogical;
57 typedef char logical1;
58 typedef char integer1;
63 /* Extern is for use with -E */
74 /*external read, write*/
83 /*internal read, write*/
113 /*rewind, backspace, endfile*/
125 ftnint *inex; /*parameters in standard's order*/
151 union Multitype { /* for multiple entry points */
162 typedef union Multitype Multitype;
164 struct Vardesc { /* for Namelist */
170 typedef struct Vardesc Vardesc;
177 typedef struct Namelist Namelist;
179 #define abs(x) ((x) >= 0 ? (x) : -(x))
180 #define dabs(x) (fabs(x))
181 #define f2cmin(a,b) ((a) <= (b) ? (a) : (b))
182 #define f2cmax(a,b) ((a) >= (b) ? (a) : (b))
183 #define dmin(a,b) (f2cmin(a,b))
184 #define dmax(a,b) (f2cmax(a,b))
185 #define bit_test(a,b) ((a) >> (b) & 1)
186 #define bit_clear(a,b) ((a) & ~((uinteger)1 << (b)))
187 #define bit_set(a,b) ((a) | ((uinteger)1 << (b)))
189 #define abort_() { sig_die("Fortran abort routine called", 1); }
190 #define c_abs(z) (cabsf(Cf(z)))
191 #define c_cos(R,Z) { pCf(R)=ccos(Cf(Z)); }
193 #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]);}
194 #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]);}
196 #define c_div(c, a, b) {pCf(c) = Cf(a)/Cf(b);}
197 #define z_div(c, a, b) {pCd(c) = Cd(a)/Cd(b);}
199 #define c_exp(R, Z) {pCf(R) = cexpf(Cf(Z));}
200 #define c_log(R, Z) {pCf(R) = clogf(Cf(Z));}
201 #define c_sin(R, Z) {pCf(R) = csinf(Cf(Z));}
202 //#define c_sqrt(R, Z) {*(R) = csqrtf(Cf(Z));}
203 #define c_sqrt(R, Z) {pCf(R) = csqrtf(Cf(Z));}
204 #define d_abs(x) (fabs(*(x)))
205 #define d_acos(x) (acos(*(x)))
206 #define d_asin(x) (asin(*(x)))
207 #define d_atan(x) (atan(*(x)))
208 #define d_atn2(x, y) (atan2(*(x),*(y)))
209 #define d_cnjg(R, Z) { pCd(R) = conj(Cd(Z)); }
210 #define r_cnjg(R, Z) { pCf(R) = conjf(Cf(Z)); }
211 #define d_cos(x) (cos(*(x)))
212 #define d_cosh(x) (cosh(*(x)))
213 #define d_dim(__a, __b) ( *(__a) > *(__b) ? *(__a) - *(__b) : 0.0 )
214 #define d_exp(x) (exp(*(x)))
215 #define d_imag(z) (cimag(Cd(z)))
216 #define r_imag(z) (cimagf(Cf(z)))
217 #define d_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
218 #define r_int(__x) (*(__x)>0 ? floor(*(__x)) : -floor(- *(__x)))
219 #define d_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
220 #define r_lg10(x) ( 0.43429448190325182765 * log(*(x)) )
221 #define d_log(x) (log(*(x)))
222 #define d_mod(x, y) (fmod(*(x), *(y)))
223 #define u_nint(__x) ((__x)>=0 ? floor((__x) + .5) : -floor(.5 - (__x)))
224 #define d_nint(x) u_nint(*(x))
225 #define u_sign(__a,__b) ((__b) >= 0 ? ((__a) >= 0 ? (__a) : -(__a)) : -((__a) >= 0 ? (__a) : -(__a)))
226 #define d_sign(a,b) u_sign(*(a),*(b))
227 #define r_sign(a,b) u_sign(*(a),*(b))
228 #define d_sin(x) (sin(*(x)))
229 #define d_sinh(x) (sinh(*(x)))
230 #define d_sqrt(x) (sqrt(*(x)))
231 #define d_tan(x) (tan(*(x)))
232 #define d_tanh(x) (tanh(*(x)))
233 #define i_abs(x) abs(*(x))
234 #define i_dnnt(x) ((integer)u_nint(*(x)))
235 #define i_len(s, n) (n)
236 #define i_nint(x) ((integer)u_nint(*(x)))
237 #define i_sign(a,b) ((integer)u_sign((integer)*(a),(integer)*(b)))
238 #define pow_dd(ap, bp) ( pow(*(ap), *(bp)))
239 #define pow_si(B,E) spow_ui(*(B),*(E))
240 #define pow_ri(B,E) spow_ui(*(B),*(E))
241 #define pow_di(B,E) dpow_ui(*(B),*(E))
242 #define pow_zi(p, a, b) {pCd(p) = zpow_ui(Cd(a), *(b));}
243 #define pow_ci(p, a, b) {pCf(p) = cpow_ui(Cf(a), *(b));}
244 #define pow_zz(R,A,B) {pCd(R) = cpow(Cd(A),*(B));}
245 #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++ = ' '; }
246 #define s_cmp(a,b,c,d) ((integer)strncmp((a),(b),f2cmin((c),(d))))
247 #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]; }
248 #define sig_die(s, kill) { exit(1); }
249 #define s_stop(s, n) {exit(0);}
250 static char junk[] = "\n@(#)LIBF77 VERSION 19990503\n";
251 #define z_abs(z) (cabs(Cd(z)))
252 #define z_exp(R, Z) {pCd(R) = cexp(Cd(Z));}
253 #define z_sqrt(R, Z) {pCd(R) = csqrt(Cd(Z));}
254 #define myexit_() break;
255 #define mycycle() continue;
256 #define myceiling(w) {ceil(w)}
257 #define myhuge(w) {HUGE_VAL}
258 //#define mymaxloc_(w,s,e,n) {if (sizeof(*(w)) == sizeof(double)) dmaxloc_((w),*(s),*(e),n); else dmaxloc_((w),*(s),*(e),n);}
259 #define mymaxloc(w,s,e,n) {dmaxloc_(w,*(s),*(e),n)}
261 /* procedure parameter types for -A and -C++ */
263 #define F2C_proc_par_types 1
265 typedef logical (*L_fp)(...);
267 typedef logical (*L_fp)();
270 static float spow_ui(float x, integer n) {
271 float pow=1.0; unsigned long int u;
273 if(n < 0) n = -n, x = 1/x;
282 static double dpow_ui(double x, integer n) {
283 double pow=1.0; unsigned long int u;
285 if(n < 0) n = -n, x = 1/x;
295 static _Fcomplex cpow_ui(complex x, integer n) {
296 complex pow={1.0,0.0}; unsigned long int u;
298 if(n < 0) n = -n, x.r = 1/x.r, x.i=1/x.i;
300 if(u & 01) pow.r *= x.r, pow.i *= x.i;
301 if(u >>= 1) x.r *= x.r, x.i *= x.i;
305 _Fcomplex p={pow.r, pow.i};
309 static _Complex float cpow_ui(_Complex float x, integer n) {
310 _Complex float pow=1.0; unsigned long int u;
312 if(n < 0) n = -n, x = 1/x;
323 static _Dcomplex zpow_ui(_Dcomplex x, integer n) {
324 _Dcomplex pow={1.0,0.0}; unsigned long int u;
326 if(n < 0) n = -n, x._Val[0] = 1/x._Val[0], x._Val[1] =1/x._Val[1];
328 if(u & 01) pow._Val[0] *= x._Val[0], pow._Val[1] *= x._Val[1];
329 if(u >>= 1) x._Val[0] *= x._Val[0], x._Val[1] *= x._Val[1];
333 _Dcomplex p = {pow._Val[0], pow._Val[1]};
337 static _Complex double zpow_ui(_Complex double x, integer n) {
338 _Complex double pow=1.0; unsigned long int u;
340 if(n < 0) n = -n, x = 1/x;
350 static integer pow_ii(integer x, integer n) {
351 integer pow; unsigned long int u;
353 if (n == 0 || x == 1) pow = 1;
354 else if (x != -1) pow = x == 0 ? 1/x : 0;
357 if ((n > 0) || !(n == 0 || x == 1 || x != -1)) {
367 static integer dmaxloc_(double *w, integer s, integer e, integer *n)
369 double m; integer i, mi;
370 for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
371 if (w[i-1]>m) mi=i ,m=w[i-1];
374 static integer smaxloc_(float *w, integer s, integer e, integer *n)
376 float m; integer i, mi;
377 for(m=w[s-1], mi=s, i=s+1; i<=e; i++)
378 if (w[i-1]>m) mi=i ,m=w[i-1];
381 static inline void cdotc_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
382 integer n = *n_, incx = *incx_, incy = *incy_, i;
384 _Fcomplex zdotc = {0.0, 0.0};
385 if (incx == 1 && incy == 1) {
386 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
387 zdotc._Val[0] += conjf(Cf(&x[i]))._Val[0] * Cf(&y[i])._Val[0];
388 zdotc._Val[1] += conjf(Cf(&x[i]))._Val[1] * Cf(&y[i])._Val[1];
391 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
392 zdotc._Val[0] += conjf(Cf(&x[i*incx]))._Val[0] * Cf(&y[i*incy])._Val[0];
393 zdotc._Val[1] += conjf(Cf(&x[i*incx]))._Val[1] * Cf(&y[i*incy])._Val[1];
399 _Complex float zdotc = 0.0;
400 if (incx == 1 && incy == 1) {
401 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
402 zdotc += conjf(Cf(&x[i])) * Cf(&y[i]);
405 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
406 zdotc += conjf(Cf(&x[i*incx])) * Cf(&y[i*incy]);
412 static inline void zdotc_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
413 integer n = *n_, incx = *incx_, incy = *incy_, i;
415 _Dcomplex zdotc = {0.0, 0.0};
416 if (incx == 1 && incy == 1) {
417 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
418 zdotc._Val[0] += conj(Cd(&x[i]))._Val[0] * Cd(&y[i])._Val[0];
419 zdotc._Val[1] += conj(Cd(&x[i]))._Val[1] * Cd(&y[i])._Val[1];
422 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
423 zdotc._Val[0] += conj(Cd(&x[i*incx]))._Val[0] * Cd(&y[i*incy])._Val[0];
424 zdotc._Val[1] += conj(Cd(&x[i*incx]))._Val[1] * Cd(&y[i*incy])._Val[1];
430 _Complex double zdotc = 0.0;
431 if (incx == 1 && incy == 1) {
432 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
433 zdotc += conj(Cd(&x[i])) * Cd(&y[i]);
436 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
437 zdotc += conj(Cd(&x[i*incx])) * Cd(&y[i*incy]);
443 static inline void cdotu_(complex *z, integer *n_, complex *x, integer *incx_, complex *y, integer *incy_) {
444 integer n = *n_, incx = *incx_, incy = *incy_, i;
446 _Fcomplex zdotc = {0.0, 0.0};
447 if (incx == 1 && incy == 1) {
448 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
449 zdotc._Val[0] += Cf(&x[i])._Val[0] * Cf(&y[i])._Val[0];
450 zdotc._Val[1] += Cf(&x[i])._Val[1] * Cf(&y[i])._Val[1];
453 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
454 zdotc._Val[0] += Cf(&x[i*incx])._Val[0] * Cf(&y[i*incy])._Val[0];
455 zdotc._Val[1] += Cf(&x[i*incx])._Val[1] * Cf(&y[i*incy])._Val[1];
461 _Complex float zdotc = 0.0;
462 if (incx == 1 && incy == 1) {
463 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
464 zdotc += Cf(&x[i]) * Cf(&y[i]);
467 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
468 zdotc += Cf(&x[i*incx]) * Cf(&y[i*incy]);
474 static inline void zdotu_(doublecomplex *z, integer *n_, doublecomplex *x, integer *incx_, doublecomplex *y, integer *incy_) {
475 integer n = *n_, incx = *incx_, incy = *incy_, i;
477 _Dcomplex zdotc = {0.0, 0.0};
478 if (incx == 1 && incy == 1) {
479 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
480 zdotc._Val[0] += Cd(&x[i])._Val[0] * Cd(&y[i])._Val[0];
481 zdotc._Val[1] += Cd(&x[i])._Val[1] * Cd(&y[i])._Val[1];
484 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
485 zdotc._Val[0] += Cd(&x[i*incx])._Val[0] * Cd(&y[i*incy])._Val[0];
486 zdotc._Val[1] += Cd(&x[i*incx])._Val[1] * Cd(&y[i*incy])._Val[1];
492 _Complex double zdotc = 0.0;
493 if (incx == 1 && incy == 1) {
494 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
495 zdotc += Cd(&x[i]) * Cd(&y[i]);
498 for (i=0;i<n;i++) { /* zdotc = zdotc + dconjg(x(i))* y(i) */
499 zdotc += Cd(&x[i*incx]) * Cd(&y[i*incy]);
505 /* -- translated by f2c (version 20000121).
506 You must link the resulting object file with the libraries:
507 -lf2c -lm (in that order)
513 /* Table of constant values */
515 static doublereal c_b12 = 0.;
516 static doublereal c_b13 = 1.;
517 static integer c__1 = 1;
518 static integer c__3 = 3;
520 /* > \brief \b DHGEQZ */
522 /* =========== DOCUMENTATION =========== */
524 /* Online html documentation available at */
525 /* http://www.netlib.org/lapack/explore-html/ */
528 /* > Download DHGEQZ + dependencies */
529 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dhgeqz.
532 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dhgeqz.
535 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dhgeqz.
543 /* SUBROUTINE DHGEQZ( JOB, COMPQ, COMPZ, N, ILO, IHI, H, LDH, T, LDT, */
544 /* ALPHAR, ALPHAI, BETA, Q, LDQ, Z, LDZ, WORK, */
547 /* CHARACTER COMPQ, COMPZ, JOB */
548 /* INTEGER IHI, ILO, INFO, LDH, LDQ, LDT, LDZ, LWORK, N */
549 /* DOUBLE PRECISION ALPHAI( * ), ALPHAR( * ), BETA( * ), */
550 /* $ H( LDH, * ), Q( LDQ, * ), T( LDT, * ), */
551 /* $ WORK( * ), Z( LDZ, * ) */
554 /* > \par Purpose: */
559 /* > DHGEQZ computes the eigenvalues of a real matrix pair (H,T), */
560 /* > where H is an upper Hessenberg matrix and T is upper triangular, */
561 /* > using the double-shift QZ method. */
562 /* > Matrix pairs of this type are produced by the reduction to */
563 /* > generalized upper Hessenberg form of a real matrix pair (A,B): */
565 /* > A = Q1*H*Z1**T, B = Q1*T*Z1**T, */
567 /* > as computed by DGGHRD. */
569 /* > If JOB='S', then the Hessenberg-triangular pair (H,T) is */
570 /* > also reduced to generalized Schur form, */
572 /* > H = Q*S*Z**T, T = Q*P*Z**T, */
574 /* > where Q and Z are orthogonal matrices, P is an upper triangular */
575 /* > matrix, and S is a quasi-triangular matrix with 1-by-1 and 2-by-2 */
576 /* > diagonal blocks. */
578 /* > The 1-by-1 blocks correspond to real eigenvalues of the matrix pair */
579 /* > (H,T) and the 2-by-2 blocks correspond to complex conjugate pairs of */
582 /* > Additionally, the 2-by-2 upper triangular diagonal blocks of P */
583 /* > corresponding to 2-by-2 blocks of S are reduced to positive diagonal */
584 /* > form, i.e., if S(j+1,j) is non-zero, then P(j+1,j) = P(j,j+1) = 0, */
585 /* > P(j,j) > 0, and P(j+1,j+1) > 0. */
587 /* > Optionally, the orthogonal matrix Q from the generalized Schur */
588 /* > factorization may be postmultiplied into an input matrix Q1, and the */
589 /* > orthogonal matrix Z may be postmultiplied into an input matrix Z1. */
590 /* > If Q1 and Z1 are the orthogonal matrices from DGGHRD that reduced */
591 /* > the matrix pair (A,B) to generalized upper Hessenberg form, then the */
592 /* > output matrices Q1*Q and Z1*Z are the orthogonal factors from the */
593 /* > generalized Schur factorization of (A,B): */
595 /* > A = (Q1*Q)*S*(Z1*Z)**T, B = (Q1*Q)*P*(Z1*Z)**T. */
597 /* > To avoid overflow, eigenvalues of the matrix pair (H,T) (equivalently, */
598 /* > of (A,B)) are computed as a pair of values (alpha,beta), where alpha is */
599 /* > complex and beta real. */
600 /* > If beta is nonzero, lambda = alpha / beta is an eigenvalue of the */
601 /* > generalized nonsymmetric eigenvalue problem (GNEP) */
602 /* > A*x = lambda*B*x */
603 /* > and if alpha is nonzero, mu = beta / alpha is an eigenvalue of the */
604 /* > alternate form of the GNEP */
605 /* > mu*A*y = B*y. */
606 /* > Real eigenvalues can be read directly from the generalized Schur */
608 /* > alpha = S(i,i), beta = P(i,i). */
610 /* > Ref: C.B. Moler & G.W. Stewart, "An Algorithm for Generalized Matrix */
611 /* > Eigenvalue Problems", SIAM J. Numer. Anal., 10(1973), */
612 /* > pp. 241--256. */
618 /* > \param[in] JOB */
620 /* > JOB is CHARACTER*1 */
621 /* > = 'E': Compute eigenvalues only; */
622 /* > = 'S': Compute eigenvalues and the Schur form. */
625 /* > \param[in] COMPQ */
627 /* > COMPQ is CHARACTER*1 */
628 /* > = 'N': Left Schur vectors (Q) are not computed; */
629 /* > = 'I': Q is initialized to the unit matrix and the matrix Q */
630 /* > of left Schur vectors of (H,T) is returned; */
631 /* > = 'V': Q must contain an orthogonal matrix Q1 on entry and */
632 /* > the product Q1*Q is returned. */
635 /* > \param[in] COMPZ */
637 /* > COMPZ is CHARACTER*1 */
638 /* > = 'N': Right Schur vectors (Z) are not computed; */
639 /* > = 'I': Z is initialized to the unit matrix and the matrix Z */
640 /* > of right Schur vectors of (H,T) is returned; */
641 /* > = 'V': Z must contain an orthogonal matrix Z1 on entry and */
642 /* > the product Z1*Z is returned. */
648 /* > The order of the matrices H, T, Q, and Z. N >= 0. */
651 /* > \param[in] ILO */
653 /* > ILO is INTEGER */
656 /* > \param[in] IHI */
658 /* > IHI is INTEGER */
659 /* > ILO and IHI mark the rows and columns of H which are in */
660 /* > Hessenberg form. It is assumed that A is already upper */
661 /* > triangular in rows and columns 1:ILO-1 and IHI+1:N. */
662 /* > If N > 0, 1 <= ILO <= IHI <= N; if N = 0, ILO=1 and IHI=0. */
665 /* > \param[in,out] H */
667 /* > H is DOUBLE PRECISION array, dimension (LDH, N) */
668 /* > On entry, the N-by-N upper Hessenberg matrix H. */
669 /* > On exit, if JOB = 'S', H contains the upper quasi-triangular */
670 /* > matrix S from the generalized Schur factorization. */
671 /* > If JOB = 'E', the diagonal blocks of H match those of S, but */
672 /* > the rest of H is unspecified. */
675 /* > \param[in] LDH */
677 /* > LDH is INTEGER */
678 /* > The leading dimension of the array H. LDH >= f2cmax( 1, N ). */
681 /* > \param[in,out] T */
683 /* > T is DOUBLE PRECISION array, dimension (LDT, N) */
684 /* > On entry, the N-by-N upper triangular matrix T. */
685 /* > On exit, if JOB = 'S', T contains the upper triangular */
686 /* > matrix P from the generalized Schur factorization; */
687 /* > 2-by-2 diagonal blocks of P corresponding to 2-by-2 blocks of S */
688 /* > are reduced to positive diagonal form, i.e., if H(j+1,j) is */
689 /* > non-zero, then T(j+1,j) = T(j,j+1) = 0, T(j,j) > 0, and */
690 /* > T(j+1,j+1) > 0. */
691 /* > If JOB = 'E', the diagonal blocks of T match those of P, but */
692 /* > the rest of T is unspecified. */
695 /* > \param[in] LDT */
697 /* > LDT is INTEGER */
698 /* > The leading dimension of the array T. LDT >= f2cmax( 1, N ). */
701 /* > \param[out] ALPHAR */
703 /* > ALPHAR is DOUBLE PRECISION array, dimension (N) */
704 /* > The real parts of each scalar alpha defining an eigenvalue */
708 /* > \param[out] ALPHAI */
710 /* > ALPHAI is DOUBLE PRECISION array, dimension (N) */
711 /* > The imaginary parts of each scalar alpha defining an */
712 /* > eigenvalue of GNEP. */
713 /* > If ALPHAI(j) is zero, then the j-th eigenvalue is real; if */
714 /* > positive, then the j-th and (j+1)-st eigenvalues are a */
715 /* > complex conjugate pair, with ALPHAI(j+1) = -ALPHAI(j). */
718 /* > \param[out] BETA */
720 /* > BETA is DOUBLE PRECISION array, dimension (N) */
721 /* > The scalars beta that define the eigenvalues of GNEP. */
722 /* > Together, the quantities alpha = (ALPHAR(j),ALPHAI(j)) and */
723 /* > beta = BETA(j) represent the j-th eigenvalue of the matrix */
724 /* > pair (A,B), in one of the forms lambda = alpha/beta or */
725 /* > mu = beta/alpha. Since either lambda or mu may overflow, */
726 /* > they should not, in general, be computed. */
729 /* > \param[in,out] Q */
731 /* > Q is DOUBLE PRECISION array, dimension (LDQ, N) */
732 /* > On entry, if COMPQ = 'V', the orthogonal matrix Q1 used in */
733 /* > the reduction of (A,B) to generalized Hessenberg form. */
734 /* > On exit, if COMPQ = 'I', the orthogonal matrix of left Schur */
735 /* > vectors of (H,T), and if COMPQ = 'V', the orthogonal matrix */
736 /* > of left Schur vectors of (A,B). */
737 /* > Not referenced if COMPQ = 'N'. */
740 /* > \param[in] LDQ */
742 /* > LDQ is INTEGER */
743 /* > The leading dimension of the array Q. LDQ >= 1. */
744 /* > If COMPQ='V' or 'I', then LDQ >= N. */
747 /* > \param[in,out] Z */
749 /* > Z is DOUBLE PRECISION array, dimension (LDZ, N) */
750 /* > On entry, if COMPZ = 'V', the orthogonal matrix Z1 used in */
751 /* > the reduction of (A,B) to generalized Hessenberg form. */
752 /* > On exit, if COMPZ = 'I', the orthogonal matrix of */
753 /* > right Schur vectors of (H,T), and if COMPZ = 'V', the */
754 /* > orthogonal matrix of right Schur vectors of (A,B). */
755 /* > Not referenced if COMPZ = 'N'. */
758 /* > \param[in] LDZ */
760 /* > LDZ is INTEGER */
761 /* > The leading dimension of the array Z. LDZ >= 1. */
762 /* > If COMPZ='V' or 'I', then LDZ >= N. */
765 /* > \param[out] WORK */
767 /* > WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK)) */
768 /* > On exit, if INFO >= 0, WORK(1) returns the optimal LWORK. */
771 /* > \param[in] LWORK */
773 /* > LWORK is INTEGER */
774 /* > The dimension of the array WORK. LWORK >= f2cmax(1,N). */
776 /* > If LWORK = -1, then a workspace query is assumed; the routine */
777 /* > only calculates the optimal size of the WORK array, returns */
778 /* > this value as the first entry of the WORK array, and no error */
779 /* > message related to LWORK is issued by XERBLA. */
782 /* > \param[out] INFO */
784 /* > INFO is INTEGER */
785 /* > = 0: successful exit */
786 /* > < 0: if INFO = -i, the i-th argument had an illegal value */
787 /* > = 1,...,N: the QZ iteration did not converge. (H,T) is not */
788 /* > in Schur form, but ALPHAR(i), ALPHAI(i), and */
789 /* > BETA(i), i=INFO+1,...,N should be correct. */
790 /* > = N+1,...,2*N: the shift calculation failed. (H,T) is not */
791 /* > in Schur form, but ALPHAR(i), ALPHAI(i), and */
792 /* > BETA(i), i=INFO-N+1,...,N should be correct. */
798 /* > \author Univ. of Tennessee */
799 /* > \author Univ. of California Berkeley */
800 /* > \author Univ. of Colorado Denver */
801 /* > \author NAG Ltd. */
803 /* > \date June 2016 */
805 /* > \ingroup doubleGEcomputational */
807 /* > \par Further Details: */
808 /* ===================== */
812 /* > Iteration counters: */
814 /* > JITER -- counts iterations. */
815 /* > IITER -- counts iterations run since ILAST was last */
816 /* > changed. This is therefore reset only when a 1-by-1 or */
817 /* > 2-by-2 block deflates off the bottom. */
820 /* ===================================================================== */
821 /* Subroutine */ int dhgeqz_(char *job, char *compq, char *compz, integer *n,
822 integer *ilo, integer *ihi, doublereal *h__, integer *ldh, doublereal
823 *t, integer *ldt, doublereal *alphar, doublereal *alphai, doublereal *
824 beta, doublereal *q, integer *ldq, doublereal *z__, integer *ldz,
825 doublereal *work, integer *lwork, integer *info)
827 /* System generated locals */
828 integer h_dim1, h_offset, q_dim1, q_offset, t_dim1, t_offset, z_dim1,
829 z_offset, i__1, i__2, i__3, i__4;
830 doublereal d__1, d__2, d__3, d__4;
832 /* Local variables */
833 doublereal ad11l, ad12l, ad21l, ad22l, ad32l, wabs, atol, btol, temp;
834 extern /* Subroutine */ int drot_(integer *, doublereal *, integer *,
835 doublereal *, integer *, doublereal *, doublereal *), dlag2_(
836 doublereal *, integer *, doublereal *, integer *, doublereal *,
837 doublereal *, doublereal *, doublereal *, doublereal *,
839 doublereal temp2, s1inv, c__;
841 doublereal s, v[3], scale;
842 extern logical lsame_(char *, char *);
843 integer iiter, ilast, jiter;
844 doublereal anorm, bnorm;
846 doublereal tempi, tempr, s1, s2, t1, u1, u2;
847 extern doublereal dlapy2_(doublereal *, doublereal *), dlapy3_(doublereal
848 *, doublereal *, doublereal *);
849 extern /* Subroutine */ int dlasv2_(doublereal *, doublereal *,
850 doublereal *, doublereal *, doublereal *, doublereal *,
851 doublereal *, doublereal *, doublereal *);
853 doublereal a11, a12, a21, a22, b11, b22, c12, c21;
855 doublereal an, bn, cl, cq, cr;
857 doublereal ascale, bscale, u12, w11;
859 doublereal cz, sl, w12, w21, w22, wi;
860 extern doublereal dlamch_(char *);
862 extern /* Subroutine */ int dlarfg_(integer *, doublereal *, doublereal *,
863 integer *, doublereal *);
865 extern doublereal dlanhs_(char *, integer *, doublereal *, integer *,
867 extern /* Subroutine */ int dlaset_(char *, integer *, integer *,
868 doublereal *, doublereal *, doublereal *, integer *);
870 extern /* Subroutine */ int dlartg_(doublereal *, doublereal *,
871 doublereal *, doublereal *, doublereal *);
873 extern /* Subroutine */ int xerbla_(char *, integer *, ftnlen);
877 integer icompq, ilastm;
882 integer icompz, ifirst;
888 doublereal a2r, b1r, b2r;
890 doublereal wr2, ad11, ad12, ad21, ad22, c11i, c22i;
892 doublereal c11r, c22r;
894 doublereal u12l, tau, sqi;
896 doublereal ulp, sqr, szi, szr;
899 /* -- LAPACK computational routine (version 3.7.0) -- */
900 /* -- LAPACK is a software package provided by Univ. of Tennessee, -- */
901 /* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
905 /* ===================================================================== */
907 /* $ SAFETY = 1.0E+0 ) */
909 /* Decode JOB, COMPQ, COMPZ */
911 /* Parameter adjustments */
913 h_offset = 1 + h_dim1 * 1;
916 t_offset = 1 + t_dim1 * 1;
922 q_offset = 1 + q_dim1 * 1;
925 z_offset = 1 + z_dim1 * 1;
930 if (lsame_(job, "E")) {
933 } else if (lsame_(job, "S")) {
940 if (lsame_(compq, "N")) {
943 } else if (lsame_(compq, "V")) {
946 } else if (lsame_(compq, "I")) {
953 if (lsame_(compz, "N")) {
956 } else if (lsame_(compz, "V")) {
959 } else if (lsame_(compz, "I")) {
966 /* Check Argument Values */
969 work[1] = (doublereal) f2cmax(1,*n);
970 lquery = *lwork == -1;
973 } else if (icompq == 0) {
975 } else if (icompz == 0) {
979 } else if (*ilo < 1) {
981 } else if (*ihi > *n || *ihi < *ilo - 1) {
983 } else if (*ldh < *n) {
985 } else if (*ldt < *n) {
987 } else if (*ldq < 1 || ilq && *ldq < *n) {
989 } else if (*ldz < 1 || ilz && *ldz < *n) {
991 } else if (*lwork < f2cmax(1,*n) && ! lquery) {
996 xerbla_("DHGEQZ", &i__1, (ftnlen)6);
1002 /* Quick return if possible */
1009 /* Initialize Q and Z */
1012 dlaset_("Full", n, n, &c_b12, &c_b13, &q[q_offset], ldq);
1015 dlaset_("Full", n, n, &c_b12, &c_b13, &z__[z_offset], ldz);
1018 /* Machine Constants */
1020 in = *ihi + 1 - *ilo;
1021 safmin = dlamch_("S");
1022 safmax = 1. / safmin;
1023 ulp = dlamch_("E") * dlamch_("B");
1024 anorm = dlanhs_("F", &in, &h__[*ilo + *ilo * h_dim1], ldh, &work[1]);
1025 bnorm = dlanhs_("F", &in, &t[*ilo + *ilo * t_dim1], ldt, &work[1]);
1027 d__1 = safmin, d__2 = ulp * anorm;
1028 atol = f2cmax(d__1,d__2);
1030 d__1 = safmin, d__2 = ulp * bnorm;
1031 btol = f2cmax(d__1,d__2);
1032 ascale = 1. / f2cmax(safmin,anorm);
1033 bscale = 1. / f2cmax(safmin,bnorm);
1035 /* Set Eigenvalues IHI+1:N */
1038 for (j = *ihi + 1; j <= i__1; ++j) {
1039 if (t[j + j * t_dim1] < 0.) {
1042 for (jr = 1; jr <= i__2; ++jr) {
1043 h__[jr + j * h_dim1] = -h__[jr + j * h_dim1];
1044 t[jr + j * t_dim1] = -t[jr + j * t_dim1];
1048 h__[j + j * h_dim1] = -h__[j + j * h_dim1];
1049 t[j + j * t_dim1] = -t[j + j * t_dim1];
1053 for (jr = 1; jr <= i__2; ++jr) {
1054 z__[jr + j * z_dim1] = -z__[jr + j * z_dim1];
1059 alphar[j] = h__[j + j * h_dim1];
1061 beta[j] = t[j + j * t_dim1];
1065 /* If IHI < ILO, skip QZ steps */
1071 /* MAIN QZ ITERATION LOOP */
1073 /* Initialize dynamic indices */
1075 /* Eigenvalues ILAST+1:N have been found. */
1076 /* Column operations modify rows IFRSTM:whatever. */
1077 /* Row operations modify columns whatever:ILASTM. */
1079 /* If only eigenvalues are being computed, then */
1080 /* IFRSTM is the row of the last splitting row above row ILAST; */
1081 /* this is always at least ILO. */
1082 /* IITER counts iterations since the last eigenvalue was found, */
1083 /* to tell when to use an extraordinary shift. */
1084 /* MAXIT is the maximum number of QZ sweeps allowed. */
1096 maxit = (*ihi - *ilo + 1) * 30;
1099 for (jiter = 1; jiter <= i__1; ++jiter) {
1101 /* Split the matrix if possible. */
1104 /* 1: H(j,j-1)=0 or j=ILO */
1107 if (ilast == *ilo) {
1109 /* Special case: j=ILAST */
1113 if ((d__1 = h__[ilast + (ilast - 1) * h_dim1], abs(d__1)) <= atol)
1115 h__[ilast + (ilast - 1) * h_dim1] = 0.;
1120 if ((d__1 = t[ilast + ilast * t_dim1], abs(d__1)) <= btol) {
1121 t[ilast + ilast * t_dim1] = 0.;
1125 /* General case: j<ILAST */
1128 for (j = ilast - 1; j >= i__2; --j) {
1130 /* Test 1: for H(j,j-1)=0 or j=ILO */
1135 if ((d__1 = h__[j + (j - 1) * h_dim1], abs(d__1)) <= atol) {
1136 h__[j + (j - 1) * h_dim1] = 0.;
1143 /* Test 2: for T(j,j)=0 */
1145 if ((d__1 = t[j + j * t_dim1], abs(d__1)) < btol) {
1146 t[j + j * t_dim1] = 0.;
1148 /* Test 1a: Check for 2 consecutive small subdiagonals in A */
1152 temp = (d__1 = h__[j + (j - 1) * h_dim1], abs(d__1));
1153 temp2 = (d__1 = h__[j + j * h_dim1], abs(d__1));
1154 tempr = f2cmax(temp,temp2);
1155 if (tempr < 1. && tempr != 0.) {
1159 if (temp * (ascale * (d__1 = h__[j + 1 + j * h_dim1], abs(
1160 d__1))) <= temp2 * (ascale * atol)) {
1165 /* If both tests pass (1 & 2), i.e., the leading diagonal */
1166 /* element of B in the block is zero, split a 1x1 block off */
1167 /* at the top. (I.e., at the J-th row/column) The leading */
1168 /* diagonal element of the remainder can also be zero, so */
1169 /* this may have to be done repeatedly. */
1171 if (ilazro || ilazr2) {
1173 for (jch = j; jch <= i__3; ++jch) {
1174 temp = h__[jch + jch * h_dim1];
1175 dlartg_(&temp, &h__[jch + 1 + jch * h_dim1], &c__, &s,
1176 &h__[jch + jch * h_dim1]);
1177 h__[jch + 1 + jch * h_dim1] = 0.;
1178 i__4 = ilastm - jch;
1179 drot_(&i__4, &h__[jch + (jch + 1) * h_dim1], ldh, &
1180 h__[jch + 1 + (jch + 1) * h_dim1], ldh, &c__,
1182 i__4 = ilastm - jch;
1183 drot_(&i__4, &t[jch + (jch + 1) * t_dim1], ldt, &t[
1184 jch + 1 + (jch + 1) * t_dim1], ldt, &c__, &s);
1186 drot_(n, &q[jch * q_dim1 + 1], &c__1, &q[(jch + 1)
1187 * q_dim1 + 1], &c__1, &c__, &s);
1190 h__[jch + (jch - 1) * h_dim1] *= c__;
1193 if ((d__1 = t[jch + 1 + (jch + 1) * t_dim1], abs(d__1)
1195 if (jch + 1 >= ilast) {
1202 t[jch + 1 + (jch + 1) * t_dim1] = 0.;
1208 /* Only test 2 passed -- chase the zero to T(ILAST,ILAST) */
1209 /* Then process as in the case T(ILAST,ILAST)=0 */
1212 for (jch = j; jch <= i__3; ++jch) {
1213 temp = t[jch + (jch + 1) * t_dim1];
1214 dlartg_(&temp, &t[jch + 1 + (jch + 1) * t_dim1], &c__,
1215 &s, &t[jch + (jch + 1) * t_dim1]);
1216 t[jch + 1 + (jch + 1) * t_dim1] = 0.;
1217 if (jch < ilastm - 1) {
1218 i__4 = ilastm - jch - 1;
1219 drot_(&i__4, &t[jch + (jch + 2) * t_dim1], ldt, &
1220 t[jch + 1 + (jch + 2) * t_dim1], ldt, &
1223 i__4 = ilastm - jch + 2;
1224 drot_(&i__4, &h__[jch + (jch - 1) * h_dim1], ldh, &
1225 h__[jch + 1 + (jch - 1) * h_dim1], ldh, &c__,
1228 drot_(n, &q[jch * q_dim1 + 1], &c__1, &q[(jch + 1)
1229 * q_dim1 + 1], &c__1, &c__, &s);
1231 temp = h__[jch + 1 + jch * h_dim1];
1232 dlartg_(&temp, &h__[jch + 1 + (jch - 1) * h_dim1], &
1233 c__, &s, &h__[jch + 1 + jch * h_dim1]);
1234 h__[jch + 1 + (jch - 1) * h_dim1] = 0.;
1235 i__4 = jch + 1 - ifrstm;
1236 drot_(&i__4, &h__[ifrstm + jch * h_dim1], &c__1, &h__[
1237 ifrstm + (jch - 1) * h_dim1], &c__1, &c__, &s)
1239 i__4 = jch - ifrstm;
1240 drot_(&i__4, &t[ifrstm + jch * t_dim1], &c__1, &t[
1241 ifrstm + (jch - 1) * t_dim1], &c__1, &c__, &s)
1244 drot_(n, &z__[jch * z_dim1 + 1], &c__1, &z__[(jch
1245 - 1) * z_dim1 + 1], &c__1, &c__, &s);
1251 } else if (ilazro) {
1253 /* Only test 1 passed -- work on J:ILAST */
1259 /* Neither test passed -- try next J */
1264 /* (Drop-through is "impossible") */
1269 /* T(ILAST,ILAST)=0 -- clear H(ILAST,ILAST-1) to split off a */
1273 temp = h__[ilast + ilast * h_dim1];
1274 dlartg_(&temp, &h__[ilast + (ilast - 1) * h_dim1], &c__, &s, &h__[
1275 ilast + ilast * h_dim1]);
1276 h__[ilast + (ilast - 1) * h_dim1] = 0.;
1277 i__2 = ilast - ifrstm;
1278 drot_(&i__2, &h__[ifrstm + ilast * h_dim1], &c__1, &h__[ifrstm + (
1279 ilast - 1) * h_dim1], &c__1, &c__, &s);
1280 i__2 = ilast - ifrstm;
1281 drot_(&i__2, &t[ifrstm + ilast * t_dim1], &c__1, &t[ifrstm + (ilast -
1282 1) * t_dim1], &c__1, &c__, &s);
1284 drot_(n, &z__[ilast * z_dim1 + 1], &c__1, &z__[(ilast - 1) *
1285 z_dim1 + 1], &c__1, &c__, &s);
1288 /* H(ILAST,ILAST-1)=0 -- Standardize B, set ALPHAR, ALPHAI, */
1292 if (t[ilast + ilast * t_dim1] < 0.) {
1295 for (j = ifrstm; j <= i__2; ++j) {
1296 h__[j + ilast * h_dim1] = -h__[j + ilast * h_dim1];
1297 t[j + ilast * t_dim1] = -t[j + ilast * t_dim1];
1301 h__[ilast + ilast * h_dim1] = -h__[ilast + ilast * h_dim1];
1302 t[ilast + ilast * t_dim1] = -t[ilast + ilast * t_dim1];
1306 for (j = 1; j <= i__2; ++j) {
1307 z__[j + ilast * z_dim1] = -z__[j + ilast * z_dim1];
1312 alphar[ilast] = h__[ilast + ilast * h_dim1];
1314 beta[ilast] = t[ilast + ilast * t_dim1];
1316 /* Go to next block -- exit if finished. */
1323 /* Reset counters */
1329 if (ifrstm > ilast) {
1337 /* This iteration only involves rows/columns IFIRST:ILAST. We */
1338 /* assume IFIRST < ILAST, and that the diagonal of B is non-zero. */
1346 /* Compute single shifts. */
1348 /* At this point, IFIRST < ILAST, and the diagonal elements of */
1349 /* T(IFIRST:ILAST,IFIRST,ILAST) are larger than BTOL (in */
1352 if (iiter / 10 * 10 == iiter) {
1354 /* Exceptional shift. Chosen for no particularly good reason. */
1355 /* (Single shift only.) */
1357 if ((doublereal) maxit * safmin * (d__1 = h__[ilast + (ilast - 1)
1358 * h_dim1], abs(d__1)) < (d__2 = t[ilast - 1 + (ilast - 1)
1359 * t_dim1], abs(d__2))) {
1360 eshift = h__[ilast + (ilast - 1) * h_dim1] / t[ilast - 1 + (
1361 ilast - 1) * t_dim1];
1363 eshift += 1. / (safmin * (doublereal) maxit);
1370 /* Shifts based on the generalized eigenvalues of the */
1371 /* bottom-right 2x2 block of A and B. The first eigenvalue */
1372 /* returned by DLAG2 is the Wilkinson shift (AEP p.512), */
1374 d__1 = safmin * 100.;
1375 dlag2_(&h__[ilast - 1 + (ilast - 1) * h_dim1], ldh, &t[ilast - 1
1376 + (ilast - 1) * t_dim1], ldt, &d__1, &s1, &s2, &wr, &wr2,
1379 if ((d__1 = wr / s1 * t[ilast + ilast * t_dim1] - h__[ilast +
1380 ilast * h_dim1], abs(d__1)) > (d__2 = wr2 / s2 * t[ilast
1381 + ilast * t_dim1] - h__[ilast + ilast * h_dim1], abs(d__2)
1392 d__3 = 1., d__4 = abs(wr), d__3 = f2cmax(d__3,d__4), d__4 = abs(wi);
1393 d__1 = s1, d__2 = safmin * f2cmax(d__3,d__4);
1394 temp = f2cmax(d__1,d__2);
1400 /* Fiddle with shift to avoid overflow */
1402 temp = f2cmin(ascale,1.) * (safmax * .5);
1409 temp = f2cmin(bscale,1.) * (safmax * .5);
1410 if (abs(wr) > temp) {
1412 d__1 = scale, d__2 = temp / abs(wr);
1413 scale = f2cmin(d__1,d__2);
1418 /* Now check for two consecutive small subdiagonals. */
1421 for (j = ilast - 1; j >= i__2; --j) {
1423 temp = (d__1 = s1 * h__[j + (j - 1) * h_dim1], abs(d__1));
1424 temp2 = (d__1 = s1 * h__[j + j * h_dim1] - wr * t[j + j * t_dim1],
1426 tempr = f2cmax(temp,temp2);
1427 if (tempr < 1. && tempr != 0.) {
1431 if ((d__1 = ascale * h__[j + 1 + j * h_dim1] * temp, abs(d__1)) <=
1432 ascale * atol * temp2) {
1441 /* Do an implicit single-shift QZ sweep. */
1445 temp = s1 * h__[istart + istart * h_dim1] - wr * t[istart + istart *
1447 temp2 = s1 * h__[istart + 1 + istart * h_dim1];
1448 dlartg_(&temp, &temp2, &c__, &s, &tempr);
1453 for (j = istart; j <= i__2; ++j) {
1455 temp = h__[j + (j - 1) * h_dim1];
1456 dlartg_(&temp, &h__[j + 1 + (j - 1) * h_dim1], &c__, &s, &h__[
1457 j + (j - 1) * h_dim1]);
1458 h__[j + 1 + (j - 1) * h_dim1] = 0.;
1462 for (jc = j; jc <= i__3; ++jc) {
1463 temp = c__ * h__[j + jc * h_dim1] + s * h__[j + 1 + jc *
1465 h__[j + 1 + jc * h_dim1] = -s * h__[j + jc * h_dim1] + c__ *
1466 h__[j + 1 + jc * h_dim1];
1467 h__[j + jc * h_dim1] = temp;
1468 temp2 = c__ * t[j + jc * t_dim1] + s * t[j + 1 + jc * t_dim1];
1469 t[j + 1 + jc * t_dim1] = -s * t[j + jc * t_dim1] + c__ * t[j
1471 t[j + jc * t_dim1] = temp2;
1476 for (jr = 1; jr <= i__3; ++jr) {
1477 temp = c__ * q[jr + j * q_dim1] + s * q[jr + (j + 1) *
1479 q[jr + (j + 1) * q_dim1] = -s * q[jr + j * q_dim1] + c__ *
1480 q[jr + (j + 1) * q_dim1];
1481 q[jr + j * q_dim1] = temp;
1486 temp = t[j + 1 + (j + 1) * t_dim1];
1487 dlartg_(&temp, &t[j + 1 + j * t_dim1], &c__, &s, &t[j + 1 + (j +
1489 t[j + 1 + j * t_dim1] = 0.;
1493 i__3 = f2cmin(i__4,ilast);
1494 for (jr = ifrstm; jr <= i__3; ++jr) {
1495 temp = c__ * h__[jr + (j + 1) * h_dim1] + s * h__[jr + j *
1497 h__[jr + j * h_dim1] = -s * h__[jr + (j + 1) * h_dim1] + c__ *
1498 h__[jr + j * h_dim1];
1499 h__[jr + (j + 1) * h_dim1] = temp;
1503 for (jr = ifrstm; jr <= i__3; ++jr) {
1504 temp = c__ * t[jr + (j + 1) * t_dim1] + s * t[jr + j * t_dim1]
1506 t[jr + j * t_dim1] = -s * t[jr + (j + 1) * t_dim1] + c__ * t[
1508 t[jr + (j + 1) * t_dim1] = temp;
1513 for (jr = 1; jr <= i__3; ++jr) {
1514 temp = c__ * z__[jr + (j + 1) * z_dim1] + s * z__[jr + j *
1516 z__[jr + j * z_dim1] = -s * z__[jr + (j + 1) * z_dim1] +
1517 c__ * z__[jr + j * z_dim1];
1518 z__[jr + (j + 1) * z_dim1] = temp;
1527 /* Use Francis double-shift */
1529 /* Note: the Francis double-shift should work with real shifts, */
1530 /* but only if the block is at least 3x3. */
1531 /* This code may break if this point is reached with */
1532 /* a 2x2 block with real eigenvalues. */
1535 if (ifirst + 1 == ilast) {
1537 /* Special case -- 2x2 block with complex eigenvectors */
1539 /* Step 1: Standardize, that is, rotate so that */
1542 /* B = ( ) with B11 non-negative. */
1545 dlasv2_(&t[ilast - 1 + (ilast - 1) * t_dim1], &t[ilast - 1 +
1546 ilast * t_dim1], &t[ilast + ilast * t_dim1], &b22, &b11, &
1556 i__2 = ilastm + 1 - ifirst;
1557 drot_(&i__2, &h__[ilast - 1 + (ilast - 1) * h_dim1], ldh, &h__[
1558 ilast + (ilast - 1) * h_dim1], ldh, &cl, &sl);
1559 i__2 = ilast + 1 - ifrstm;
1560 drot_(&i__2, &h__[ifrstm + (ilast - 1) * h_dim1], &c__1, &h__[
1561 ifrstm + ilast * h_dim1], &c__1, &cr, &sr);
1563 if (ilast < ilastm) {
1564 i__2 = ilastm - ilast;
1565 drot_(&i__2, &t[ilast - 1 + (ilast + 1) * t_dim1], ldt, &t[
1566 ilast + (ilast + 1) * t_dim1], ldt, &cl, &sl);
1568 if (ifrstm < ilast - 1) {
1569 i__2 = ifirst - ifrstm;
1570 drot_(&i__2, &t[ifrstm + (ilast - 1) * t_dim1], &c__1, &t[
1571 ifrstm + ilast * t_dim1], &c__1, &cr, &sr);
1575 drot_(n, &q[(ilast - 1) * q_dim1 + 1], &c__1, &q[ilast *
1576 q_dim1 + 1], &c__1, &cl, &sl);
1579 drot_(n, &z__[(ilast - 1) * z_dim1 + 1], &c__1, &z__[ilast *
1580 z_dim1 + 1], &c__1, &cr, &sr);
1583 t[ilast - 1 + (ilast - 1) * t_dim1] = b11;
1584 t[ilast - 1 + ilast * t_dim1] = 0.;
1585 t[ilast + (ilast - 1) * t_dim1] = 0.;
1586 t[ilast + ilast * t_dim1] = b22;
1588 /* If B22 is negative, negate column ILAST */
1592 for (j = ifrstm; j <= i__2; ++j) {
1593 h__[j + ilast * h_dim1] = -h__[j + ilast * h_dim1];
1594 t[j + ilast * t_dim1] = -t[j + ilast * t_dim1];
1600 for (j = 1; j <= i__2; ++j) {
1601 z__[j + ilast * z_dim1] = -z__[j + ilast * z_dim1];
1608 /* Step 2: Compute ALPHAR, ALPHAI, and BETA (see refs.) */
1610 /* Recompute shift */
1612 d__1 = safmin * 100.;
1613 dlag2_(&h__[ilast - 1 + (ilast - 1) * h_dim1], ldh, &t[ilast - 1
1614 + (ilast - 1) * t_dim1], ldt, &d__1, &s1, &temp, &wr, &
1617 /* If standardization has perturbed the shift onto real line, */
1618 /* do another (real single-shift) QR step. */
1625 /* Do EISPACK (QZVAL) computation of alpha and beta */
1627 a11 = h__[ilast - 1 + (ilast - 1) * h_dim1];
1628 a21 = h__[ilast + (ilast - 1) * h_dim1];
1629 a12 = h__[ilast - 1 + ilast * h_dim1];
1630 a22 = h__[ilast + ilast * h_dim1];
1632 /* Compute complex Givens rotation on right */
1633 /* (Assume some element of C = (sA - wB) > unfl ) */
1635 /* (sA - wB) ( CZ -SZ ) */
1638 c11r = s1 * a11 - wr * b11;
1642 c22r = s1 * a22 - wr * b22;
1645 if (abs(c11r) + abs(c11i) + abs(c12) > abs(c21) + abs(c22r) + abs(
1647 t1 = dlapy3_(&c12, &c11r, &c11i);
1652 cz = dlapy2_(&c22r, &c22i);
1660 t1 = dlapy2_(&cz, &c21);
1662 szr = -c21 * tempr / t1;
1663 szi = c21 * tempi / t1;
1667 /* Compute Givens rotation on left */
1673 an = abs(a11) + abs(a12) + abs(a21) + abs(a22);
1674 bn = abs(b11) + abs(b22);
1675 wabs = abs(wr) + abs(wi);
1676 if (s1 * an > wabs * bn) {
1681 a1r = cz * a11 + szr * a12;
1683 a2r = cz * a21 + szr * a22;
1685 cq = dlapy2_(&a1r, &a1i);
1693 sqr = tempr * a2r + tempi * a2i;
1694 sqi = tempi * a2r - tempr * a2i;
1697 t1 = dlapy3_(&cq, &sqr, &sqi);
1702 /* Compute diagonal elements of QBZ */
1704 tempr = sqr * szr - sqi * szi;
1705 tempi = sqr * szi + sqi * szr;
1706 b1r = cq * cz * b11 + tempr * b22;
1708 b1a = dlapy2_(&b1r, &b1i);
1709 b2r = cq * cz * b22 + tempr * b11;
1711 b2a = dlapy2_(&b2r, &b2i);
1713 /* Normalize so beta > 0, and Im( alpha1 ) > 0 */
1715 beta[ilast - 1] = b1a;
1717 alphar[ilast - 1] = wr * b1a * s1inv;
1718 alphai[ilast - 1] = wi * b1a * s1inv;
1719 alphar[ilast] = wr * b2a * s1inv;
1720 alphai[ilast] = -(wi * b2a) * s1inv;
1722 /* Step 3: Go to next block -- exit if finished. */
1729 /* Reset counters */
1735 if (ifrstm > ilast) {
1742 /* Usual case: 3x3 or larger block, using Francis implicit */
1746 /* Eigenvalue equation is w - c w + d = 0, */
1749 /* so compute 1st column of (A B ) - c A B + d */
1750 /* using the formula in QZIT (from EISPACK) */
1752 /* We assume that the block is at least 3x3 */
1754 ad11 = ascale * h__[ilast - 1 + (ilast - 1) * h_dim1] / (bscale *
1755 t[ilast - 1 + (ilast - 1) * t_dim1]);
1756 ad21 = ascale * h__[ilast + (ilast - 1) * h_dim1] / (bscale * t[
1757 ilast - 1 + (ilast - 1) * t_dim1]);
1758 ad12 = ascale * h__[ilast - 1 + ilast * h_dim1] / (bscale * t[
1759 ilast + ilast * t_dim1]);
1760 ad22 = ascale * h__[ilast + ilast * h_dim1] / (bscale * t[ilast +
1762 u12 = t[ilast - 1 + ilast * t_dim1] / t[ilast + ilast * t_dim1];
1763 ad11l = ascale * h__[ifirst + ifirst * h_dim1] / (bscale * t[
1764 ifirst + ifirst * t_dim1]);
1765 ad21l = ascale * h__[ifirst + 1 + ifirst * h_dim1] / (bscale * t[
1766 ifirst + ifirst * t_dim1]);
1767 ad12l = ascale * h__[ifirst + (ifirst + 1) * h_dim1] / (bscale *
1768 t[ifirst + 1 + (ifirst + 1) * t_dim1]);
1769 ad22l = ascale * h__[ifirst + 1 + (ifirst + 1) * h_dim1] / (
1770 bscale * t[ifirst + 1 + (ifirst + 1) * t_dim1]);
1771 ad32l = ascale * h__[ifirst + 2 + (ifirst + 1) * h_dim1] / (
1772 bscale * t[ifirst + 1 + (ifirst + 1) * t_dim1]);
1773 u12l = t[ifirst + (ifirst + 1) * t_dim1] / t[ifirst + 1 + (ifirst
1776 v[0] = (ad11 - ad11l) * (ad22 - ad11l) - ad12 * ad21 + ad21 * u12
1777 * ad11l + (ad12l - ad11l * u12l) * ad21l;
1778 v[1] = (ad22l - ad11l - ad21l * u12l - (ad11 - ad11l) - (ad22 -
1779 ad11l) + ad21 * u12) * ad21l;
1780 v[2] = ad32l * ad21l;
1784 dlarfg_(&c__3, v, &v[1], &c__1, &tau);
1790 for (j = istart; j <= i__2; ++j) {
1792 /* All but last elements: use 3x3 Householder transforms. */
1794 /* Zero (j-1)st column of A */
1797 v[0] = h__[j + (j - 1) * h_dim1];
1798 v[1] = h__[j + 1 + (j - 1) * h_dim1];
1799 v[2] = h__[j + 2 + (j - 1) * h_dim1];
1801 dlarfg_(&c__3, &h__[j + (j - 1) * h_dim1], &v[1], &c__1, &
1804 h__[j + 1 + (j - 1) * h_dim1] = 0.;
1805 h__[j + 2 + (j - 1) * h_dim1] = 0.;
1809 for (jc = j; jc <= i__3; ++jc) {
1810 temp = tau * (h__[j + jc * h_dim1] + v[1] * h__[j + 1 +
1811 jc * h_dim1] + v[2] * h__[j + 2 + jc * h_dim1]);
1812 h__[j + jc * h_dim1] -= temp;
1813 h__[j + 1 + jc * h_dim1] -= temp * v[1];
1814 h__[j + 2 + jc * h_dim1] -= temp * v[2];
1815 temp2 = tau * (t[j + jc * t_dim1] + v[1] * t[j + 1 + jc *
1816 t_dim1] + v[2] * t[j + 2 + jc * t_dim1]);
1817 t[j + jc * t_dim1] -= temp2;
1818 t[j + 1 + jc * t_dim1] -= temp2 * v[1];
1819 t[j + 2 + jc * t_dim1] -= temp2 * v[2];
1824 for (jr = 1; jr <= i__3; ++jr) {
1825 temp = tau * (q[jr + j * q_dim1] + v[1] * q[jr + (j +
1826 1) * q_dim1] + v[2] * q[jr + (j + 2) * q_dim1]
1828 q[jr + j * q_dim1] -= temp;
1829 q[jr + (j + 1) * q_dim1] -= temp * v[1];
1830 q[jr + (j + 2) * q_dim1] -= temp * v[2];
1835 /* Zero j-th column of B (see DLAGBC for details) */
1837 /* Swap rows to pivot */
1841 d__3 = (d__1 = t[j + 1 + (j + 1) * t_dim1], abs(d__1)), d__4 =
1842 (d__2 = t[j + 1 + (j + 2) * t_dim1], abs(d__2));
1843 temp = f2cmax(d__3,d__4);
1845 d__3 = (d__1 = t[j + 2 + (j + 1) * t_dim1], abs(d__1)), d__4 =
1846 (d__2 = t[j + 2 + (j + 2) * t_dim1], abs(d__2));
1847 temp2 = f2cmax(d__3,d__4);
1848 if (f2cmax(temp,temp2) < safmin) {
1853 } else if (temp >= temp2) {
1854 w11 = t[j + 1 + (j + 1) * t_dim1];
1855 w21 = t[j + 2 + (j + 1) * t_dim1];
1856 w12 = t[j + 1 + (j + 2) * t_dim1];
1857 w22 = t[j + 2 + (j + 2) * t_dim1];
1858 u1 = t[j + 1 + j * t_dim1];
1859 u2 = t[j + 2 + j * t_dim1];
1861 w21 = t[j + 1 + (j + 1) * t_dim1];
1862 w11 = t[j + 2 + (j + 1) * t_dim1];
1863 w22 = t[j + 1 + (j + 2) * t_dim1];
1864 w12 = t[j + 2 + (j + 2) * t_dim1];
1865 u2 = t[j + 1 + j * t_dim1];
1866 u1 = t[j + 2 + j * t_dim1];
1869 /* Swap columns if nec. */
1871 if (abs(w12) > abs(w11)) {
1891 if (abs(w22) < safmin) {
1897 if (abs(w22) < abs(u2)) {
1898 scale = (d__1 = w22 / u2, abs(d__1));
1900 if (abs(w11) < abs(u1)) {
1902 d__2 = scale, d__3 = (d__1 = w11 / u1, abs(d__1));
1903 scale = f2cmin(d__2,d__3);
1908 u2 = scale * u2 / w22;
1909 u1 = (scale * u1 - w12 * u2) / w11;
1918 /* Compute Householder Vector */
1920 /* Computing 2nd power */
1922 /* Computing 2nd power */
1924 /* Computing 2nd power */
1926 t1 = sqrt(d__1 * d__1 + d__2 * d__2 + d__3 * d__3);
1927 tau = scale / t1 + 1.;
1928 vs = -1. / (scale + t1);
1933 /* Apply transformations from the right. */
1937 i__3 = f2cmin(i__4,ilast);
1938 for (jr = ifrstm; jr <= i__3; ++jr) {
1939 temp = tau * (h__[jr + j * h_dim1] + v[1] * h__[jr + (j +
1940 1) * h_dim1] + v[2] * h__[jr + (j + 2) * h_dim1]);
1941 h__[jr + j * h_dim1] -= temp;
1942 h__[jr + (j + 1) * h_dim1] -= temp * v[1];
1943 h__[jr + (j + 2) * h_dim1] -= temp * v[2];
1947 for (jr = ifrstm; jr <= i__3; ++jr) {
1948 temp = tau * (t[jr + j * t_dim1] + v[1] * t[jr + (j + 1) *
1949 t_dim1] + v[2] * t[jr + (j + 2) * t_dim1]);
1950 t[jr + j * t_dim1] -= temp;
1951 t[jr + (j + 1) * t_dim1] -= temp * v[1];
1952 t[jr + (j + 2) * t_dim1] -= temp * v[2];
1957 for (jr = 1; jr <= i__3; ++jr) {
1958 temp = tau * (z__[jr + j * z_dim1] + v[1] * z__[jr + (
1959 j + 1) * z_dim1] + v[2] * z__[jr + (j + 2) *
1961 z__[jr + j * z_dim1] -= temp;
1962 z__[jr + (j + 1) * z_dim1] -= temp * v[1];
1963 z__[jr + (j + 2) * z_dim1] -= temp * v[2];
1967 t[j + 1 + j * t_dim1] = 0.;
1968 t[j + 2 + j * t_dim1] = 0.;
1972 /* Last elements: Use Givens rotations */
1974 /* Rotations from the left */
1977 temp = h__[j + (j - 1) * h_dim1];
1978 dlartg_(&temp, &h__[j + 1 + (j - 1) * h_dim1], &c__, &s, &h__[j +
1980 h__[j + 1 + (j - 1) * h_dim1] = 0.;
1983 for (jc = j; jc <= i__2; ++jc) {
1984 temp = c__ * h__[j + jc * h_dim1] + s * h__[j + 1 + jc *
1986 h__[j + 1 + jc * h_dim1] = -s * h__[j + jc * h_dim1] + c__ *
1987 h__[j + 1 + jc * h_dim1];
1988 h__[j + jc * h_dim1] = temp;
1989 temp2 = c__ * t[j + jc * t_dim1] + s * t[j + 1 + jc * t_dim1];
1990 t[j + 1 + jc * t_dim1] = -s * t[j + jc * t_dim1] + c__ * t[j
1992 t[j + jc * t_dim1] = temp2;
1997 for (jr = 1; jr <= i__2; ++jr) {
1998 temp = c__ * q[jr + j * q_dim1] + s * q[jr + (j + 1) *
2000 q[jr + (j + 1) * q_dim1] = -s * q[jr + j * q_dim1] + c__ *
2001 q[jr + (j + 1) * q_dim1];
2002 q[jr + j * q_dim1] = temp;
2007 /* Rotations from the right. */
2009 temp = t[j + 1 + (j + 1) * t_dim1];
2010 dlartg_(&temp, &t[j + 1 + j * t_dim1], &c__, &s, &t[j + 1 + (j +
2012 t[j + 1 + j * t_dim1] = 0.;
2015 for (jr = ifrstm; jr <= i__2; ++jr) {
2016 temp = c__ * h__[jr + (j + 1) * h_dim1] + s * h__[jr + j *
2018 h__[jr + j * h_dim1] = -s * h__[jr + (j + 1) * h_dim1] + c__ *
2019 h__[jr + j * h_dim1];
2020 h__[jr + (j + 1) * h_dim1] = temp;
2024 for (jr = ifrstm; jr <= i__2; ++jr) {
2025 temp = c__ * t[jr + (j + 1) * t_dim1] + s * t[jr + j * t_dim1]
2027 t[jr + j * t_dim1] = -s * t[jr + (j + 1) * t_dim1] + c__ * t[
2029 t[jr + (j + 1) * t_dim1] = temp;
2034 for (jr = 1; jr <= i__2; ++jr) {
2035 temp = c__ * z__[jr + (j + 1) * z_dim1] + s * z__[jr + j *
2037 z__[jr + j * z_dim1] = -s * z__[jr + (j + 1) * z_dim1] +
2038 c__ * z__[jr + j * z_dim1];
2039 z__[jr + (j + 1) * z_dim1] = temp;
2044 /* End of Double-Shift code */
2050 /* End of iteration loop */
2057 /* Drop-through = non-convergence */
2062 /* Successful completion of all QZ steps */
2066 /* Set Eigenvalues 1:ILO-1 */
2069 for (j = 1; j <= i__1; ++j) {
2070 if (t[j + j * t_dim1] < 0.) {
2073 for (jr = 1; jr <= i__2; ++jr) {
2074 h__[jr + j * h_dim1] = -h__[jr + j * h_dim1];
2075 t[jr + j * t_dim1] = -t[jr + j * t_dim1];
2079 h__[j + j * h_dim1] = -h__[j + j * h_dim1];
2080 t[j + j * t_dim1] = -t[j + j * t_dim1];
2084 for (jr = 1; jr <= i__2; ++jr) {
2085 z__[jr + j * z_dim1] = -z__[jr + j * z_dim1];
2090 alphar[j] = h__[j + j * h_dim1];
2092 beta[j] = t[j + j * t_dim1];
2096 /* Normal Termination */
2100 /* Exit (other than argument error) -- return optimal workspace size */
2103 work[1] = (doublereal) (*n);