5 #include "tnt_math_utils.h"
36 for (
int j = 0; j < n; j++) {
42 for (
int i = n-1; i > 0; i--) {
48 for (
int k = 0; k < i; k++) {
49 scale = scale + abs(d[k]);
53 for (
int j = 0; j < i; j++) {
62 for (
int k = 0; k < i; k++) {
74 for (
int j = 0; j < i; j++) {
80 for (
int j = 0; j < i; j++) {
83 g = e[j] + V[j][j] * f;
84 for (
int k = j+1; k <= i-1; k++) {
91 for (
int j = 0; j < i; j++) {
95 float hh = f / (h + h);
96 for (
int j = 0; j < i; j++) {
99 for (
int j = 0; j < i; j++) {
102 for (
int k = j; k <= i-1; k++) {
103 V[k][j] -= (f * e[k] + g * d[k]);
114 for (
int i = 0; i < n-1; i++) {
119 for (
int k = 0; k <= i; k++) {
120 d[k] = V[k][i+1] / h;
122 for (
int j = 0; j <= i; j++) {
124 for (
int k = 0; k <= i; k++) {
125 g += V[k][i+1] * V[k][j];
127 for (
int k = 0; k <= i; k++) {
132 for (
int k = 0; k <= i; k++) {
136 for (
int j = 0; j < n; j++) {
153 for (
int i = 1; i < n; i++) {
160 float eps = pow(2.0,-52.0);
161 for (
int l = 0; l < n; l++) {
165 tst1 = max(tst1,abs(d[l]) + abs(e[l]));
170 if (abs(e[m]) <= eps*tst1) {
188 float p = (d[l+1] - g) / (2.0 * e[l]);
189 float r = sqrt( p*p + 1 );
193 d[l] = e[l] / (p + r);
194 d[l+1] = e[l] * (p + r);
197 for (
int i = l+2; i < n; i++) {
211 for (
int i = m-1; i >= l; i--) {
217 r = sqrt( p*p + e[i]*e[i] );
221 p = c * d[i] - s * g;
222 d[i+1] = h + s * (c * g + s * d[i]);
226 for (
int k = 0; k < n; k++) {
228 V[k][i+1] = s * V[k][i] + c * h;
229 V[k][i] = c * V[k][i] - s * h;
232 p = -s * s2 * c3 * el1 * e[l] / dl1;
238 }
while (abs(e[l]) > eps*tst1);
246 for (
int i = 0; i < n-1; i++) {
249 for (
int j = i+1; j < n; j++) {
258 for (
int j = 0; j < n; j++) {
279 for (
int m = low+1; m <= high-1; m++) {
284 for (
int i = m; i <= high; i++) {
285 scale = scale + abs(H[i][m-1]);
292 for (
int i = high; i >= m; i--) {
293 ort[i] = H[i][m-1]/scale;
294 h += ort[i] * ort[i];
306 for (
int j = m; j < n; j++) {
308 for (
int i = high; i >= m; i--) {
312 for (
int i = m; i <= high; i++) {
317 for (
int i = 0; i <= high; i++) {
319 for (
int j = high; j >= m; j--) {
323 for (
int j = m; j <= high; j++) {
327 ort[m] = scale*ort[m];
334 for (
int i = 0; i < n; i++) {
335 for (
int j = 0; j < n; j++) {
336 V[i][j] = (i == j ? 1.0 : 0.0);
340 for (
int m = high-1; m >= low+1; m--) {
341 if (H[m][m-1] != 0.0) {
342 for (
int i = m+1; i <= high; i++) {
345 for (
int j = m; j <= high; j++) {
347 for (
int i = m; i <= high; i++) {
348 g += ort[i] * V[i][j];
351 g = (g / ort[m]) / H[m][m-1];
352 for (
int i = m; i <= high; i++) {
353 V[i][j] += g * ort[i];
363 float cdivr{}, cdivi{};
364 void cdiv(
float xr,
float xi,
float yr,
float yi) {
366 if (abs(yr) > abs(yi)) {
369 cdivr = (xr + r*xi)/d;
370 cdivi = (xi - r*xr)/d;
374 cdivr = (r*xr + xi)/d;
375 cdivi = (r*xi - xr)/d;
395 float eps = pow(2.0,-52.0);
397 float p=0,q=0,r=0,s=0,z=0,t,w,x,y;
402 for (
int i = 0; i < nn; i++) {
403 if ((i < low) || (i > high)) {
407 for (
int j = max(i-1,0); j < nn; j++) {
408 norm = norm + abs(H[i][j]);
421 s = abs(H[l-1][l-1]) + abs(H[l][l]);
425 if (abs(H[l][l-1]) < eps * s) {
435 H[n][n] = H[n][n] + exshift;
443 }
else if (l == n-1) {
444 w = H[n][n-1] * H[n-1][n];
445 p = (H[n-1][n-1] - H[n][n]) / 2.0;
448 H[n][n] = H[n][n] + exshift;
449 H[n-1][n-1] = H[n-1][n-1] + exshift;
471 r = sqrt(p * p+q * q);
477 for (
int j = n-1; j < nn; j++) {
479 H[n-1][j] = q * z + p * H[n][j];
480 H[n][j] = q * H[n][j] - p * z;
485 for (
int i = 0; i <= n; i++) {
487 H[i][n-1] = q * z + p * H[i][n];
488 H[i][n] = q * H[i][n] - p * z;
493 for (
int i = low; i <= high; i++) {
495 V[i][n-1] = q * z + p * V[i][n];
496 V[i][n] = q * V[i][n] - p * z;
521 w = H[n][n-1] * H[n-1][n];
528 for (
int i = low; i <= n; i++) {
531 s = abs(H[n][n-1]) + abs(H[n-1][n-2]);
546 s = x - w / ((y - x) / 2.0 + s);
547 for (
int i = low; i <= n; i++) {
566 p = (r * s - w) / H[m+1][m] + H[m][m+1];
567 q = H[m+1][m+1] - z - r - s;
569 s = abs(p) + abs(q) + abs(r);
578 if (abs(H[m][m-1]) * (abs(q) + abs(r)) <
579 eps * (abs(p) * (abs(H[m-1][m-1]) + abs(z) +
580 abs(H[m+1][m+1])))) {
586 for (
int i = m+2; i <= n; i++) {
595 for (
int k = m; k <= n-1; k++) {
596 int notlast = (k != n-1);
600 r = (notlast ? H[k+2][k-1] : 0.0);
601 x = abs(p) + abs(q) + abs(r);
611 s = sqrt(p * p + q * q + r * r);
619 H[k][k-1] = -H[k][k-1];
630 for (
int j = k; j < nn; j++) {
631 p = H[k][j] + q * H[k+1][j];
633 p = p + r * H[k+2][j];
634 H[k+2][j] = H[k+2][j] - p * z;
636 H[k][j] = H[k][j] - p * x;
637 H[k+1][j] = H[k+1][j] - p * y;
642 for (
int i = 0; i <= min(n,k+3); i++) {
643 p = x * H[i][k] + y * H[i][k+1];
645 p = p + z * H[i][k+2];
646 H[i][k+2] = H[i][k+2] - p * r;
648 H[i][k] = H[i][k] - p;
649 H[i][k+1] = H[i][k+1] - p * q;
654 for (
int i = low; i <= high; i++) {
655 p = x * V[i][k] + y * V[i][k+1];
657 p = p + z * V[i][k+2];
658 V[i][k+2] = V[i][k+2] - p * r;
660 V[i][k] = V[i][k] - p;
661 V[i][k+1] = V[i][k+1] - p * q;
674 for (n = nn-1; n >= 0; n--) {
683 for (
int i = n-1; i >= 0; i--) {
686 for (
int j = l; j <= n; j++) {
687 r = r + H[i][j] * H[j][n];
698 H[i][n] = -r / (eps * norm);
706 q = (d[i] - p) * (d[i] - p) + e[i] * e[i];
707 t = (x * s - z * r) / q;
709 if (abs(x) > abs(z)) {
710 H[i+1][n] = (-r - w * t) / x;
712 H[i+1][n] = (-s - y * t) / z;
719 if ((eps * t) * t > 1) {
720 for (
int j = i; j <= n; j++) {
721 H[j][n] = H[j][n] / t;
734 if (abs(H[n][n-1]) > abs(H[n-1][n])) {
735 H[n-1][n-1] = q / H[n][n-1];
736 H[n-1][n] = -(H[n][n] - p) / H[n][n-1];
738 cdiv(0.0,-H[n-1][n],H[n-1][n-1]-p,q);
744 for (
int i = n-2; i >= 0; i--) {
748 for (
int j = l; j <= n; j++) {
749 ra = ra + H[i][j] * H[j][n-1];
750 sa = sa + H[i][j] * H[j][n];
770 vr = (d[i] - p) * (d[i] - p) + e[i] * e[i] - q * q;
771 vi = (d[i] - p) * 2.0 * q;
772 if ((vr == 0.0) && (vi == 0.0)) {
773 vr = eps * norm * (abs(w) + abs(q) +
774 abs(x) + abs(y) + abs(z));
776 cdiv(x*r-z*ra+q*sa,x*s-z*sa-q*ra,vr,vi);
779 if (abs(x) > (abs(z) + abs(q))) {
780 H[i+1][n-1] = (-ra - w * H[i][n-1] + q * H[i][n]) / x;
781 H[i+1][n] = (-sa - w * H[i][n] - q * H[i][n-1]) / x;
783 cdiv(-r-y*H[i][n-1],-s-y*H[i][n],z,q);
791 t = max(abs(H[i][n-1]),abs(H[i][n]));
792 if ((eps * t) * t > 1) {
793 for (
int j = i; j <= n; j++) {
794 H[j][n-1] = H[j][n-1] / t;
795 H[j][n] = H[j][n] / t;
805 for (
int i = 0; i < nn; i++) {
808 if (i < low || i > high) {
809 for (
int j = i; j < nn; j++) {
817 for (
int j = nn-1; j >= low; j--) {
818 for (
int i = low; i <= high; i++) {
820 for (
int k = low; k <= min(j,high); k++) {
821 z = z + V[i][k] * H[k][j];
841 for (
int j = 0; (j < n) && issymmetric; j++) {
844 for (
int i = 0; (i < n) && issymmetric; i++) {
845 issymmetric = (A[i][j] == A[j][i]);
851 for (
int i = 0; i < n; i++) {
852 for (
int j = 0; j < n; j++) {
866 for (
int j = 0; j < n; j++) {
867 for (
int i = 0; i < n; i++) {
886 for (
int i=0;i<3;i++)
888 for (
int j=0;j<3;j++)
889 { V_[i][j] = V[i][j];}}
898 for (
int i=0;i<3;i++)
909 for (
int i=0;i<3;i++)
949 for (
int i = 0; i < n; i++) {
950 for (
int j = 0; j < n; j++) {
956 }
else if (e[i] < 0) {
Eigenvalue(float A[3][3])
Check for symmetry, then construct the eigenvalue decomposition.
Definition: jama_eig.h:837
void getRealEigenvalues(float d_[3])
Return the real parts of the eigenvalues.
Definition: jama_eig.h:897
void getImagEigenvalues(float e_[3])
Return the imaginary parts of the eigenvalues in parameter e_.
Definition: jama_eig.h:908
void getV(float V_[3][3])
Return the eigenvector matrix.
Definition: jama_eig.h:885
void getD(float D[3][3])
Computes the block diagonal eigenvalue matrix.
Definition: jama_eig.h:948
Definition: jama_eig.h:16