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 integer c__1 = 1;
516 static integer c__0 = 0;
517 static integer c_n1 = -1;
519 /* > \brief <b> SGEEVX computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE mat
522 /* =========== DOCUMENTATION =========== */
524 /* Online html documentation available at */
525 /* http://www.netlib.org/lapack/explore-html/ */
528 /* > Download SGEEVX + dependencies */
529 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgeevx.
532 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgeevx.
535 /* > <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgeevx.
543 /* SUBROUTINE SGEEVX( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, WR, WI, */
544 /* VL, LDVL, VR, LDVR, ILO, IHI, SCALE, ABNRM, */
545 /* RCONDE, RCONDV, WORK, LWORK, IWORK, INFO ) */
547 /* CHARACTER BALANC, JOBVL, JOBVR, SENSE */
548 /* INTEGER IHI, ILO, INFO, LDA, LDVL, LDVR, LWORK, N */
550 /* INTEGER IWORK( * ) */
551 /* REAL A( LDA, * ), RCONDE( * ), RCONDV( * ), */
552 /* $ SCALE( * ), VL( LDVL, * ), VR( LDVR, * ), */
553 /* $ WI( * ), WORK( * ), WR( * ) */
556 /* > \par Purpose: */
561 /* > SGEEVX computes for an N-by-N real nonsymmetric matrix A, the */
562 /* > eigenvalues and, optionally, the left and/or right eigenvectors. */
564 /* > Optionally also, it computes a balancing transformation to improve */
565 /* > the conditioning of the eigenvalues and eigenvectors (ILO, IHI, */
566 /* > SCALE, and ABNRM), reciprocal condition numbers for the eigenvalues */
567 /* > (RCONDE), and reciprocal condition numbers for the right */
568 /* > eigenvectors (RCONDV). */
570 /* > The right eigenvector v(j) of A satisfies */
571 /* > A * v(j) = lambda(j) * v(j) */
572 /* > where lambda(j) is its eigenvalue. */
573 /* > The left eigenvector u(j) of A satisfies */
574 /* > u(j)**H * A = lambda(j) * u(j)**H */
575 /* > where u(j)**H denotes the conjugate-transpose of u(j). */
577 /* > The computed eigenvectors are normalized to have Euclidean norm */
578 /* > equal to 1 and largest component real. */
580 /* > Balancing a matrix means permuting the rows and columns to make it */
581 /* > more nearly upper triangular, and applying a diagonal similarity */
582 /* > transformation D * A * D**(-1), where D is a diagonal matrix, to */
583 /* > make its rows and columns closer in norm and the condition numbers */
584 /* > of its eigenvalues and eigenvectors smaller. The computed */
585 /* > reciprocal condition numbers correspond to the balanced matrix. */
586 /* > Permuting rows and columns will not change the condition numbers */
587 /* > (in exact arithmetic) but diagonal scaling will. For further */
588 /* > explanation of balancing, see section 4.10.2 of the LAPACK */
589 /* > Users' Guide. */
595 /* > \param[in] BALANC */
597 /* > BALANC is CHARACTER*1 */
598 /* > Indicates how the input matrix should be diagonally scaled */
599 /* > and/or permuted to improve the conditioning of its */
601 /* > = 'N': Do not diagonally scale or permute; */
602 /* > = 'P': Perform permutations to make the matrix more nearly */
603 /* > upper triangular. Do not diagonally scale; */
604 /* > = 'S': Diagonally scale the matrix, i.e. replace A by */
605 /* > D*A*D**(-1), where D is a diagonal matrix chosen */
606 /* > to make the rows and columns of A more equal in */
607 /* > norm. Do not permute; */
608 /* > = 'B': Both diagonally scale and permute A. */
610 /* > Computed reciprocal condition numbers will be for the matrix */
611 /* > after balancing and/or permuting. Permuting does not change */
612 /* > condition numbers (in exact arithmetic), but balancing does. */
615 /* > \param[in] JOBVL */
617 /* > JOBVL is CHARACTER*1 */
618 /* > = 'N': left eigenvectors of A are not computed; */
619 /* > = 'V': left eigenvectors of A are computed. */
620 /* > If SENSE = 'E' or 'B', JOBVL must = 'V'. */
623 /* > \param[in] JOBVR */
625 /* > JOBVR is CHARACTER*1 */
626 /* > = 'N': right eigenvectors of A are not computed; */
627 /* > = 'V': right eigenvectors of A are computed. */
628 /* > If SENSE = 'E' or 'B', JOBVR must = 'V'. */
631 /* > \param[in] SENSE */
633 /* > SENSE is CHARACTER*1 */
634 /* > Determines which reciprocal condition numbers are computed. */
635 /* > = 'N': None are computed; */
636 /* > = 'E': Computed for eigenvalues only; */
637 /* > = 'V': Computed for right eigenvectors only; */
638 /* > = 'B': Computed for eigenvalues and right eigenvectors. */
640 /* > If SENSE = 'E' or 'B', both left and right eigenvectors */
641 /* > must also be computed (JOBVL = 'V' and JOBVR = 'V'). */
647 /* > The order of the matrix A. N >= 0. */
650 /* > \param[in,out] A */
652 /* > A is REAL array, dimension (LDA,N) */
653 /* > On entry, the N-by-N matrix A. */
654 /* > On exit, A has been overwritten. If JOBVL = 'V' or */
655 /* > JOBVR = 'V', A contains the real Schur form of the balanced */
656 /* > version of the input matrix A. */
659 /* > \param[in] LDA */
661 /* > LDA is INTEGER */
662 /* > The leading dimension of the array A. LDA >= f2cmax(1,N). */
665 /* > \param[out] WR */
667 /* > WR is REAL array, dimension (N) */
670 /* > \param[out] WI */
672 /* > WI is REAL array, dimension (N) */
673 /* > WR and WI contain the real and imaginary parts, */
674 /* > respectively, of the computed eigenvalues. Complex */
675 /* > conjugate pairs of eigenvalues will appear consecutively */
676 /* > with the eigenvalue having the positive imaginary part */
680 /* > \param[out] VL */
682 /* > VL is REAL array, dimension (LDVL,N) */
683 /* > If JOBVL = 'V', the left eigenvectors u(j) are stored one */
684 /* > after another in the columns of VL, in the same order */
685 /* > as their eigenvalues. */
686 /* > If JOBVL = 'N', VL is not referenced. */
687 /* > If the j-th eigenvalue is real, then u(j) = VL(:,j), */
688 /* > the j-th column of VL. */
689 /* > If the j-th and (j+1)-st eigenvalues form a complex */
690 /* > conjugate pair, then u(j) = VL(:,j) + i*VL(:,j+1) and */
691 /* > u(j+1) = VL(:,j) - i*VL(:,j+1). */
694 /* > \param[in] LDVL */
696 /* > LDVL is INTEGER */
697 /* > The leading dimension of the array VL. LDVL >= 1; if */
698 /* > JOBVL = 'V', LDVL >= N. */
701 /* > \param[out] VR */
703 /* > VR is REAL array, dimension (LDVR,N) */
704 /* > If JOBVR = 'V', the right eigenvectors v(j) are stored one */
705 /* > after another in the columns of VR, in the same order */
706 /* > as their eigenvalues. */
707 /* > If JOBVR = 'N', VR is not referenced. */
708 /* > If the j-th eigenvalue is real, then v(j) = VR(:,j), */
709 /* > the j-th column of VR. */
710 /* > If the j-th and (j+1)-st eigenvalues form a complex */
711 /* > conjugate pair, then v(j) = VR(:,j) + i*VR(:,j+1) and */
712 /* > v(j+1) = VR(:,j) - i*VR(:,j+1). */
715 /* > \param[in] LDVR */
717 /* > LDVR is INTEGER */
718 /* > The leading dimension of the array VR. LDVR >= 1, and if */
719 /* > JOBVR = 'V', LDVR >= N. */
722 /* > \param[out] ILO */
724 /* > ILO is INTEGER */
727 /* > \param[out] IHI */
729 /* > IHI is INTEGER */
730 /* > ILO and IHI are integer values determined when A was */
731 /* > balanced. The balanced A(i,j) = 0 if I > J and */
732 /* > J = 1,...,ILO-1 or I = IHI+1,...,N. */
735 /* > \param[out] SCALE */
737 /* > SCALE is REAL array, dimension (N) */
738 /* > Details of the permutations and scaling factors applied */
739 /* > when balancing A. If P(j) is the index of the row and column */
740 /* > interchanged with row and column j, and D(j) is the scaling */
741 /* > factor applied to row and column j, then */
742 /* > SCALE(J) = P(J), for J = 1,...,ILO-1 */
743 /* > = D(J), for J = ILO,...,IHI */
744 /* > = P(J) for J = IHI+1,...,N. */
745 /* > The order in which the interchanges are made is N to IHI+1, */
746 /* > then 1 to ILO-1. */
749 /* > \param[out] ABNRM */
751 /* > ABNRM is REAL */
752 /* > The one-norm of the balanced matrix (the maximum */
753 /* > of the sum of absolute values of elements of any column). */
756 /* > \param[out] RCONDE */
758 /* > RCONDE is REAL array, dimension (N) */
759 /* > RCONDE(j) is the reciprocal condition number of the j-th */
763 /* > \param[out] RCONDV */
765 /* > RCONDV is REAL array, dimension (N) */
766 /* > RCONDV(j) is the reciprocal condition number of the j-th */
767 /* > right eigenvector. */
770 /* > \param[out] WORK */
772 /* > WORK is REAL array, dimension (MAX(1,LWORK)) */
773 /* > On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
776 /* > \param[in] LWORK */
778 /* > LWORK is INTEGER */
779 /* > The dimension of the array WORK. If SENSE = 'N' or 'E', */
780 /* > LWORK >= f2cmax(1,2*N), and if JOBVL = 'V' or JOBVR = 'V', */
781 /* > LWORK >= 3*N. If SENSE = 'V' or 'B', LWORK >= N*(N+6). */
782 /* > For good performance, LWORK must generally be larger. */
784 /* > If LWORK = -1, then a workspace query is assumed; the routine */
785 /* > only calculates the optimal size of the WORK array, returns */
786 /* > this value as the first entry of the WORK array, and no error */
787 /* > message related to LWORK is issued by XERBLA. */
790 /* > \param[out] IWORK */
792 /* > IWORK is INTEGER array, dimension (2*N-2) */
793 /* > If SENSE = 'N' or 'E', not referenced. */
796 /* > \param[out] INFO */
798 /* > INFO is INTEGER */
799 /* > = 0: successful exit */
800 /* > < 0: if INFO = -i, the i-th argument had an illegal value. */
801 /* > > 0: if INFO = i, the QR algorithm failed to compute all the */
802 /* > eigenvalues, and no eigenvectors or condition numbers */
803 /* > have been computed; elements 1:ILO-1 and i+1:N of WR */
804 /* > and WI contain eigenvalues which have converged. */
810 /* > \author Univ. of Tennessee */
811 /* > \author Univ. of California Berkeley */
812 /* > \author Univ. of Colorado Denver */
813 /* > \author NAG Ltd. */
815 /* > \date June 2016 */
817 /* @generated from dgeevx.f, fortran d -> s, Tue Apr 19 01:47:44 2016 */
819 /* > \ingroup realGEeigen */
821 /* ===================================================================== */
822 /* Subroutine */ int sgeevx_(char *balanc, char *jobvl, char *jobvr, char *
823 sense, integer *n, real *a, integer *lda, real *wr, real *wi, real *
824 vl, integer *ldvl, real *vr, integer *ldvr, integer *ilo, integer *
825 ihi, real *scale, real *abnrm, real *rconde, real *rcondv, real *work,
826 integer *lwork, integer *iwork, integer *info)
828 /* System generated locals */
829 integer a_dim1, a_offset, vl_dim1, vl_offset, vr_dim1, vr_offset, i__1,
833 /* Local variables */
836 integer ierr, itau, iwrk, nout;
837 extern /* Subroutine */ int srot_(integer *, real *, integer *, real *,
838 integer *, real *, real *);
839 extern real snrm2_(integer *, real *, integer *);
843 extern logical lsame_(char *, char *);
844 extern /* Subroutine */ int sscal_(integer *, real *, real *, integer *);
845 extern real slapy2_(real *, real *);
847 extern /* Subroutine */ int slabad_(real *, real *);
850 extern /* Subroutine */ int sgebak_(char *, char *, integer *, integer *,
851 integer *, real *, integer *, real *, integer *, integer *), sgebal_(char *, integer *, real *, integer *,
852 integer *, integer *, real *, integer *);
854 extern real slamch_(char *), slange_(char *, integer *, integer *,
855 real *, integer *, real *);
856 extern /* Subroutine */ int sgehrd_(integer *, integer *, integer *, real
857 *, integer *, real *, real *, integer *, integer *), xerbla_(char
858 *, integer *, ftnlen);
859 extern integer ilaenv_(integer *, char *, char *, integer *, integer *,
860 integer *, integer *, ftnlen, ftnlen);
863 extern /* Subroutine */ int slascl_(char *, integer *, integer *, real *,
864 real *, integer *, integer *, real *, integer *, integer *);
865 extern integer isamax_(integer *, real *, integer *);
866 extern /* Subroutine */ int slacpy_(char *, integer *, integer *, real *,
867 integer *, real *, integer *), slartg_(real *, real *,
868 real *, real *, real *), sorghr_(integer *, integer *, integer *,
869 real *, integer *, real *, real *, integer *, integer *), shseqr_(
870 char *, char *, integer *, integer *, integer *, real *, integer *
871 , real *, real *, real *, integer *, real *, integer *, integer *);
872 integer minwrk, maxwrk;
873 extern /* Subroutine */ int strsna_(char *, char *, logical *, integer *,
874 real *, integer *, real *, integer *, real *, integer *, real *,
875 real *, integer *, integer *, real *, integer *, integer *,
877 logical wantvl, wntsnb;
881 logical lquery, wantvr, wntsnn, wntsnv;
882 extern /* Subroutine */ int strevc3_(char *, char *, logical *, integer *,
883 real *, integer *, real *, integer *, real *, integer *, integer
884 *, integer *, real *, integer *, integer *);
886 real scl, dum[1], eps;
887 integer lwork_trevc__;
890 /* -- LAPACK driver routine (version 3.7.1) -- */
891 /* -- LAPACK is a software package provided by Univ. of Tennessee, -- */
892 /* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..-- */
896 /* ===================================================================== */
899 /* Test the input arguments */
901 /* Parameter adjustments */
903 a_offset = 1 + a_dim1 * 1;
908 vl_offset = 1 + vl_dim1 * 1;
911 vr_offset = 1 + vr_dim1 * 1;
921 lquery = *lwork == -1;
922 wantvl = lsame_(jobvl, "V");
923 wantvr = lsame_(jobvr, "V");
924 wntsnn = lsame_(sense, "N");
925 wntsne = lsame_(sense, "E");
926 wntsnv = lsame_(sense, "V");
927 wntsnb = lsame_(sense, "B");
928 if (! (lsame_(balanc, "N") || lsame_(balanc, "S") || lsame_(balanc, "P")
929 || lsame_(balanc, "B"))) {
931 } else if (! wantvl && ! lsame_(jobvl, "N")) {
933 } else if (! wantvr && ! lsame_(jobvr, "N")) {
935 } else if (! (wntsnn || wntsne || wntsnb || wntsnv) || (wntsne || wntsnb)
936 && ! (wantvl && wantvr)) {
940 } else if (*lda < f2cmax(1,*n)) {
942 } else if (*ldvl < 1 || wantvl && *ldvl < *n) {
944 } else if (*ldvr < 1 || wantvr && *ldvr < *n) {
948 /* Compute workspace */
949 /* (Note: Comments in the code beginning "Workspace:" describe the */
950 /* minimal amount of workspace needed at that point in the code, */
951 /* as well as the preferred amount for good performance. */
952 /* NB refers to the optimal block size for the immediately */
953 /* following subroutine, as returned by ILAENV. */
954 /* HSWORK refers to the workspace preferred by SHSEQR, as */
955 /* calculated below. HSWORK is computed assuming ILO=1 and IHI=N, */
956 /* the worst case.) */
963 maxwrk = *n + *n * ilaenv_(&c__1, "SGEHRD", " ", n, &c__1, n, &
964 c__0, (ftnlen)6, (ftnlen)1);
967 strevc3_("L", "B", select, n, &a[a_offset], lda, &vl[
968 vl_offset], ldvl, &vr[vr_offset], ldvr, n, &nout, &
969 work[1], &c_n1, &ierr);
970 lwork_trevc__ = (integer) work[1];
972 i__1 = maxwrk, i__2 = *n + lwork_trevc__;
973 maxwrk = f2cmax(i__1,i__2);
974 shseqr_("S", "V", n, &c__1, n, &a[a_offset], lda, &wr[1], &wi[
975 1], &vl[vl_offset], ldvl, &work[1], &c_n1, info);
977 strevc3_("R", "B", select, n, &a[a_offset], lda, &vl[
978 vl_offset], ldvl, &vr[vr_offset], ldvr, n, &nout, &
979 work[1], &c_n1, &ierr);
980 lwork_trevc__ = (integer) work[1];
982 i__1 = maxwrk, i__2 = *n + lwork_trevc__;
983 maxwrk = f2cmax(i__1,i__2);
984 shseqr_("S", "V", n, &c__1, n, &a[a_offset], lda, &wr[1], &wi[
985 1], &vr[vr_offset], ldvr, &work[1], &c_n1, info);
988 shseqr_("E", "N", n, &c__1, n, &a[a_offset], lda, &wr[1],
989 &wi[1], &vr[vr_offset], ldvr, &work[1], &c_n1,
992 shseqr_("S", "N", n, &c__1, n, &a[a_offset], lda, &wr[1],
993 &wi[1], &vr[vr_offset], ldvr, &work[1], &c_n1,
997 hswork = (integer) work[1];
999 if (! wantvl && ! wantvr) {
1003 i__1 = minwrk, i__2 = *n * *n + *n * 6;
1004 minwrk = f2cmax(i__1,i__2);
1006 maxwrk = f2cmax(maxwrk,hswork);
1009 i__1 = maxwrk, i__2 = *n * *n + *n * 6;
1010 maxwrk = f2cmax(i__1,i__2);
1014 if (! wntsnn && ! wntsne) {
1016 i__1 = minwrk, i__2 = *n * *n + *n * 6;
1017 minwrk = f2cmax(i__1,i__2);
1019 maxwrk = f2cmax(maxwrk,hswork);
1021 i__1 = maxwrk, i__2 = *n + (*n - 1) * ilaenv_(&c__1, "SORGHR",
1022 " ", n, &c__1, n, &c_n1, (ftnlen)6, (ftnlen)1);
1023 maxwrk = f2cmax(i__1,i__2);
1024 if (! wntsnn && ! wntsne) {
1026 i__1 = maxwrk, i__2 = *n * *n + *n * 6;
1027 maxwrk = f2cmax(i__1,i__2);
1030 i__1 = maxwrk, i__2 = *n * 3;
1031 maxwrk = f2cmax(i__1,i__2);
1033 maxwrk = f2cmax(maxwrk,minwrk);
1035 work[1] = (real) maxwrk;
1037 if (*lwork < minwrk && ! lquery) {
1044 xerbla_("SGEEVX", &i__1, (ftnlen)6);
1046 } else if (lquery) {
1050 /* Quick return if possible */
1056 /* Get machine constants */
1059 smlnum = slamch_("S");
1060 bignum = 1.f / smlnum;
1061 slabad_(&smlnum, &bignum);
1062 smlnum = sqrt(smlnum) / eps;
1063 bignum = 1.f / smlnum;
1065 /* Scale A if f2cmax element outside range [SMLNUM,BIGNUM] */
1068 anrm = slange_("M", n, n, &a[a_offset], lda, dum);
1070 if (anrm > 0.f && anrm < smlnum) {
1073 } else if (anrm > bignum) {
1078 slascl_("G", &c__0, &c__0, &anrm, &cscale, n, n, &a[a_offset], lda, &
1082 /* Balance the matrix and compute ABNRM */
1084 sgebal_(balanc, n, &a[a_offset], lda, ilo, ihi, &scale[1], &ierr);
1085 *abnrm = slange_("1", n, n, &a[a_offset], lda, dum);
1088 slascl_("G", &c__0, &c__0, &cscale, &anrm, &c__1, &c__1, dum, &c__1, &
1093 /* Reduce to upper Hessenberg form */
1094 /* (Workspace: need 2*N, prefer N+N*NB) */
1098 i__1 = *lwork - iwrk + 1;
1099 sgehrd_(n, ilo, ihi, &a[a_offset], lda, &work[itau], &work[iwrk], &i__1, &
1104 /* Want left eigenvectors */
1105 /* Copy Householder vectors to VL */
1107 *(unsigned char *)side = 'L';
1108 slacpy_("L", n, n, &a[a_offset], lda, &vl[vl_offset], ldvl)
1111 /* Generate orthogonal matrix in VL */
1112 /* (Workspace: need 2*N-1, prefer N+(N-1)*NB) */
1114 i__1 = *lwork - iwrk + 1;
1115 sorghr_(n, ilo, ihi, &vl[vl_offset], ldvl, &work[itau], &work[iwrk], &
1118 /* Perform QR iteration, accumulating Schur vectors in VL */
1119 /* (Workspace: need 1, prefer HSWORK (see comments) ) */
1122 i__1 = *lwork - iwrk + 1;
1123 shseqr_("S", "V", n, ilo, ihi, &a[a_offset], lda, &wr[1], &wi[1], &vl[
1124 vl_offset], ldvl, &work[iwrk], &i__1, info);
1128 /* Want left and right eigenvectors */
1129 /* Copy Schur vectors to VR */
1131 *(unsigned char *)side = 'B';
1132 slacpy_("F", n, n, &vl[vl_offset], ldvl, &vr[vr_offset], ldvr);
1135 } else if (wantvr) {
1137 /* Want right eigenvectors */
1138 /* Copy Householder vectors to VR */
1140 *(unsigned char *)side = 'R';
1141 slacpy_("L", n, n, &a[a_offset], lda, &vr[vr_offset], ldvr)
1144 /* Generate orthogonal matrix in VR */
1145 /* (Workspace: need 2*N-1, prefer N+(N-1)*NB) */
1147 i__1 = *lwork - iwrk + 1;
1148 sorghr_(n, ilo, ihi, &vr[vr_offset], ldvr, &work[itau], &work[iwrk], &
1151 /* Perform QR iteration, accumulating Schur vectors in VR */
1152 /* (Workspace: need 1, prefer HSWORK (see comments) ) */
1155 i__1 = *lwork - iwrk + 1;
1156 shseqr_("S", "V", n, ilo, ihi, &a[a_offset], lda, &wr[1], &wi[1], &vr[
1157 vr_offset], ldvr, &work[iwrk], &i__1, info);
1161 /* Compute eigenvalues only */
1162 /* If condition numbers desired, compute Schur form */
1165 *(unsigned char *)job = 'E';
1167 *(unsigned char *)job = 'S';
1170 /* (Workspace: need 1, prefer HSWORK (see comments) ) */
1173 i__1 = *lwork - iwrk + 1;
1174 shseqr_(job, "N", n, ilo, ihi, &a[a_offset], lda, &wr[1], &wi[1], &vr[
1175 vr_offset], ldvr, &work[iwrk], &i__1, info);
1178 /* If INFO .NE. 0 from SHSEQR, then quit */
1184 if (wantvl || wantvr) {
1186 /* Compute left and/or right eigenvectors */
1187 /* (Workspace: need 3*N, prefer N + 2*N*NB) */
1189 i__1 = *lwork - iwrk + 1;
1190 strevc3_(side, "B", select, n, &a[a_offset], lda, &vl[vl_offset],
1191 ldvl, &vr[vr_offset], ldvr, n, &nout, &work[iwrk], &i__1, &
1195 /* Compute condition numbers if desired */
1196 /* (Workspace: need N*N+6*N unless SENSE = 'E') */
1199 strsna_(sense, "A", select, n, &a[a_offset], lda, &vl[vl_offset],
1200 ldvl, &vr[vr_offset], ldvr, &rconde[1], &rcondv[1], n, &nout,
1201 &work[iwrk], n, &iwork[1], &icond);
1206 /* Undo balancing of left eigenvectors */
1208 sgebak_(balanc, "L", n, ilo, ihi, &scale[1], n, &vl[vl_offset], ldvl,
1211 /* Normalize left eigenvectors and make largest component real */
1214 for (i__ = 1; i__ <= i__1; ++i__) {
1215 if (wi[i__] == 0.f) {
1216 scl = 1.f / snrm2_(n, &vl[i__ * vl_dim1 + 1], &c__1);
1217 sscal_(n, &scl, &vl[i__ * vl_dim1 + 1], &c__1);
1218 } else if (wi[i__] > 0.f) {
1219 r__1 = snrm2_(n, &vl[i__ * vl_dim1 + 1], &c__1);
1220 r__2 = snrm2_(n, &vl[(i__ + 1) * vl_dim1 + 1], &c__1);
1221 scl = 1.f / slapy2_(&r__1, &r__2);
1222 sscal_(n, &scl, &vl[i__ * vl_dim1 + 1], &c__1);
1223 sscal_(n, &scl, &vl[(i__ + 1) * vl_dim1 + 1], &c__1);
1225 for (k = 1; k <= i__2; ++k) {
1226 /* Computing 2nd power */
1227 r__1 = vl[k + i__ * vl_dim1];
1228 /* Computing 2nd power */
1229 r__2 = vl[k + (i__ + 1) * vl_dim1];
1230 work[k] = r__1 * r__1 + r__2 * r__2;
1233 k = isamax_(n, &work[1], &c__1);
1234 slartg_(&vl[k + i__ * vl_dim1], &vl[k + (i__ + 1) * vl_dim1],
1236 srot_(n, &vl[i__ * vl_dim1 + 1], &c__1, &vl[(i__ + 1) *
1237 vl_dim1 + 1], &c__1, &cs, &sn);
1238 vl[k + (i__ + 1) * vl_dim1] = 0.f;
1246 /* Undo balancing of right eigenvectors */
1248 sgebak_(balanc, "R", n, ilo, ihi, &scale[1], n, &vr[vr_offset], ldvr,
1251 /* Normalize right eigenvectors and make largest component real */
1254 for (i__ = 1; i__ <= i__1; ++i__) {
1255 if (wi[i__] == 0.f) {
1256 scl = 1.f / snrm2_(n, &vr[i__ * vr_dim1 + 1], &c__1);
1257 sscal_(n, &scl, &vr[i__ * vr_dim1 + 1], &c__1);
1258 } else if (wi[i__] > 0.f) {
1259 r__1 = snrm2_(n, &vr[i__ * vr_dim1 + 1], &c__1);
1260 r__2 = snrm2_(n, &vr[(i__ + 1) * vr_dim1 + 1], &c__1);
1261 scl = 1.f / slapy2_(&r__1, &r__2);
1262 sscal_(n, &scl, &vr[i__ * vr_dim1 + 1], &c__1);
1263 sscal_(n, &scl, &vr[(i__ + 1) * vr_dim1 + 1], &c__1);
1265 for (k = 1; k <= i__2; ++k) {
1266 /* Computing 2nd power */
1267 r__1 = vr[k + i__ * vr_dim1];
1268 /* Computing 2nd power */
1269 r__2 = vr[k + (i__ + 1) * vr_dim1];
1270 work[k] = r__1 * r__1 + r__2 * r__2;
1273 k = isamax_(n, &work[1], &c__1);
1274 slartg_(&vr[k + i__ * vr_dim1], &vr[k + (i__ + 1) * vr_dim1],
1276 srot_(n, &vr[i__ * vr_dim1 + 1], &c__1, &vr[(i__ + 1) *
1277 vr_dim1 + 1], &c__1, &cs, &sn);
1278 vr[k + (i__ + 1) * vr_dim1] = 0.f;
1284 /* Undo scaling if necessary */
1291 i__2 = f2cmax(i__3,1);
1292 slascl_("G", &c__0, &c__0, &cscale, &anrm, &i__1, &c__1, &wr[*info +
1297 i__2 = f2cmax(i__3,1);
1298 slascl_("G", &c__0, &c__0, &cscale, &anrm, &i__1, &c__1, &wi[*info +
1301 if ((wntsnv || wntsnb) && icond == 0) {
1302 slascl_("G", &c__0, &c__0, &cscale, &anrm, n, &c__1, &rcondv[
1307 slascl_("G", &c__0, &c__0, &cscale, &anrm, &i__1, &c__1, &wr[1],
1310 slascl_("G", &c__0, &c__0, &cscale, &anrm, &i__1, &c__1, &wi[1],
1315 work[1] = (real) maxwrk;