00001 #ifndef JAMA_EIG_H
00002 #define JAMA_EIG_H
00003
00004
00005 #include "tnt_array1d.h"
00006 #include "tnt_array2d.h"
00007 #include "tnt_math_utils.h"
00008
00009 #include <algorithm>
00010
00011
00012 #include <cmath>
00013
00014
00015 using namespace TNT;
00016 using namespace std;
00017
00018 namespace JAMA
00019 {
00020
00071 template <class Real>
00072 class Eigenvalue
00073 {
00074
00075
00077 int n;
00078
00079 int issymmetric;
00080
00083 TNT::Array1D<Real> d;
00084 TNT::Array1D<Real> e;
00085
00087 TNT::Array2D<Real> V;
00088
00092 TNT::Array2D<Real> H;
00093
00094
00098 TNT::Array1D<Real> ort;
00099
00100
00101
00102
00103 void tred2() {
00104
00105
00106
00107
00108
00109
00110 for (int j = 0; j < n; j++) {
00111 d[j] = V[n-1][j];
00112 }
00113
00114
00115
00116 for (int i = n-1; i > 0; i--) {
00117
00118
00119
00120 Real scale = 0.0;
00121 Real h = 0.0;
00122 for (int k = 0; k < i; k++) {
00123 scale = scale + abs(d[k]);
00124 }
00125 if (scale == 0.0) {
00126 e[i] = d[i-1];
00127 for (int j = 0; j < i; j++) {
00128 d[j] = V[i-1][j];
00129 V[i][j] = 0.0;
00130 V[j][i] = 0.0;
00131 }
00132 } else {
00133
00134
00135
00136 for (int k = 0; k < i; k++) {
00137 d[k] /= scale;
00138 h += d[k] * d[k];
00139 }
00140 Real f = d[i-1];
00141 Real g = sqrt(h);
00142 if (f > 0) {
00143 g = -g;
00144 }
00145 e[i] = scale * g;
00146 h = h - f * g;
00147 d[i-1] = f - g;
00148 for (int j = 0; j < i; j++) {
00149 e[j] = 0.0;
00150 }
00151
00152
00153
00154 for (int j = 0; j < i; j++) {
00155 f = d[j];
00156 V[j][i] = f;
00157 g = e[j] + V[j][j] * f;
00158 for (int k = j+1; k <= i-1; k++) {
00159 g += V[k][j] * d[k];
00160 e[k] += V[k][j] * f;
00161 }
00162 e[j] = g;
00163 }
00164 f = 0.0;
00165 for (int j = 0; j < i; j++) {
00166 e[j] /= h;
00167 f += e[j] * d[j];
00168 }
00169 Real hh = f / (h + h);
00170 for (int j = 0; j < i; j++) {
00171 e[j] -= hh * d[j];
00172 }
00173 for (int j = 0; j < i; j++) {
00174 f = d[j];
00175 g = e[j];
00176 for (int k = j; k <= i-1; k++) {
00177 V[k][j] -= (f * e[k] + g * d[k]);
00178 }
00179 d[j] = V[i-1][j];
00180 V[i][j] = 0.0;
00181 }
00182 }
00183 d[i] = h;
00184 }
00185
00186
00187
00188 for (int i = 0; i < n-1; i++) {
00189 V[n-1][i] = V[i][i];
00190 V[i][i] = 1.0;
00191 Real h = d[i+1];
00192 if (h != 0.0) {
00193 for (int k = 0; k <= i; k++) {
00194 d[k] = V[k][i+1] / h;
00195 }
00196 for (int j = 0; j <= i; j++) {
00197 Real g = 0.0;
00198 for (int k = 0; k <= i; k++) {
00199 g += V[k][i+1] * V[k][j];
00200 }
00201 for (int k = 0; k <= i; k++) {
00202 V[k][j] -= g * d[k];
00203 }
00204 }
00205 }
00206 for (int k = 0; k <= i; k++) {
00207 V[k][i+1] = 0.0;
00208 }
00209 }
00210 for (int j = 0; j < n; j++) {
00211 d[j] = V[n-1][j];
00212 V[n-1][j] = 0.0;
00213 }
00214 V[n-1][n-1] = 1.0;
00215 e[0] = 0.0;
00216 }
00217
00218
00219
00220 void tql2 () {
00221
00222
00223
00224
00225
00226
00227 for (int i = 1; i < n; i++) {
00228 e[i-1] = e[i];
00229 }
00230 e[n-1] = 0.0;
00231
00232 Real f = 0.0;
00233 Real tst1 = 0.0;
00234 Real eps = pow(2.0,-52.0);
00235 for (int l = 0; l < n; l++) {
00236
00237
00238
00239 tst1 = max(tst1,abs(d[l]) + abs(e[l]));
00240 int m = l;
00241
00242
00243 while (m < n) {
00244 if (abs(e[m]) <= eps*tst1) {
00245 break;
00246 }
00247 m++;
00248 }
00249
00250
00251
00252
00253
00254 if (m > l) {
00255 int iter = 0;
00256 do {
00257 iter = iter + 1;
00258
00259
00260
00261 Real g = d[l];
00262 Real p = (d[l+1] - g) / (2.0 * e[l]);
00263 Real r = hypot(p,1.0);
00264 if (p < 0) {
00265 r = -r;
00266 }
00267 d[l] = e[l] / (p + r);
00268 d[l+1] = e[l] * (p + r);
00269 Real dl1 = d[l+1];
00270 Real h = g - d[l];
00271 for (int i = l+2; i < n; i++) {
00272 d[i] -= h;
00273 }
00274 f = f + h;
00275
00276
00277
00278 p = d[m];
00279 Real c = 1.0;
00280 Real c2 = c;
00281 Real c3 = c;
00282 Real el1 = e[l+1];
00283 Real s = 0.0;
00284 Real s2 = 0.0;
00285 for (int i = m-1; i >= l; i--) {
00286 c3 = c2;
00287 c2 = c;
00288 s2 = s;
00289 g = c * e[i];
00290 h = c * p;
00291 r = hypot(p,e[i]);
00292 e[i+1] = s * r;
00293 s = e[i] / r;
00294 c = p / r;
00295 p = c * d[i] - s * g;
00296 d[i+1] = h + s * (c * g + s * d[i]);
00297
00298
00299
00300 for (int k = 0; k < n; k++) {
00301 h = V[k][i+1];
00302 V[k][i+1] = s * V[k][i] + c * h;
00303 V[k][i] = c * V[k][i] - s * h;
00304 }
00305 }
00306 p = -s * s2 * c3 * el1 * e[l] / dl1;
00307 e[l] = s * p;
00308 d[l] = c * p;
00309
00310
00311
00312 } while (abs(e[l]) > eps*tst1);
00313 }
00314 d[l] = d[l] + f;
00315 e[l] = 0.0;
00316 }
00317
00318
00319
00320 for (int i = 0; i < n-1; i++) {
00321 int k = i;
00322 Real p = d[i];
00323 for (int j = i+1; j < n; j++) {
00324 if (d[j] < p) {
00325 k = j;
00326 p = d[j];
00327 }
00328 }
00329 if (k != i) {
00330 d[k] = d[i];
00331 d[i] = p;
00332 for (int j = 0; j < n; j++) {
00333 p = V[j][i];
00334 V[j][i] = V[j][k];
00335 V[j][k] = p;
00336 }
00337 }
00338 }
00339 }
00340
00341
00342
00343 void orthes () {
00344
00345
00346
00347
00348
00349
00350 int low = 0;
00351 int high = n-1;
00352
00353 for (int m = low+1; m <= high-1; m++) {
00354
00355
00356
00357 Real scale = 0.0;
00358 for (int i = m; i <= high; i++) {
00359 scale = scale + abs(H[i][m-1]);
00360 }
00361 if (scale != 0.0) {
00362
00363
00364
00365 Real h = 0.0;
00366 for (int i = high; i >= m; i--) {
00367 ort[i] = H[i][m-1]/scale;
00368 h += ort[i] * ort[i];
00369 }
00370 Real g = sqrt(h);
00371 if (ort[m] > 0) {
00372 g = -g;
00373 }
00374 h = h - ort[m] * g;
00375 ort[m] = ort[m] - g;
00376
00377
00378
00379
00380 for (int j = m; j < n; j++) {
00381 Real f = 0.0;
00382 for (int i = high; i >= m; i--) {
00383 f += ort[i]*H[i][j];
00384 }
00385 f = f/h;
00386 for (int i = m; i <= high; i++) {
00387 H[i][j] -= f*ort[i];
00388 }
00389 }
00390
00391 for (int i = 0; i <= high; i++) {
00392 Real f = 0.0;
00393 for (int j = high; j >= m; j--) {
00394 f += ort[j]*H[i][j];
00395 }
00396 f = f/h;
00397 for (int j = m; j <= high; j++) {
00398 H[i][j] -= f*ort[j];
00399 }
00400 }
00401 ort[m] = scale*ort[m];
00402 H[m][m-1] = scale*g;
00403 }
00404 }
00405
00406
00407
00408 for (int i = 0; i < n; i++) {
00409 for (int j = 0; j < n; j++) {
00410 V[i][j] = (i == j ? 1.0 : 0.0);
00411 }
00412 }
00413
00414 for (int m = high-1; m >= low+1; m--) {
00415 if (H[m][m-1] != 0.0) {
00416 for (int i = m+1; i <= high; i++) {
00417 ort[i] = H[i][m-1];
00418 }
00419 for (int j = m; j <= high; j++) {
00420 Real g = 0.0;
00421 for (int i = m; i <= high; i++) {
00422 g += ort[i] * V[i][j];
00423 }
00424
00425 g = (g / ort[m]) / H[m][m-1];
00426 for (int i = m; i <= high; i++) {
00427 V[i][j] += g * ort[i];
00428 }
00429 }
00430 }
00431 }
00432 }
00433
00434
00435
00436
00437 Real cdivr, cdivi;
00438 void cdiv(Real xr, Real xi, Real yr, Real yi) {
00439 Real r,d;
00440 if (abs(yr) > abs(yi)) {
00441 r = yi/yr;
00442 d = yr + r*yi;
00443 cdivr = (xr + r*xi)/d;
00444 cdivi = (xi - r*xr)/d;
00445 } else {
00446 r = yr/yi;
00447 d = yi + r*yr;
00448 cdivr = (r*xr + xi)/d;
00449 cdivi = (r*xi - xr)/d;
00450 }
00451 }
00452
00453
00454
00455
00456 void hqr2 () {
00457
00458
00459
00460
00461
00462
00463
00464
00465 int nn = this->n;
00466 int n = nn-1;
00467 int low = 0;
00468 int high = nn-1;
00469 Real eps = pow(2.0,-52.0);
00470 Real exshift = 0.0;
00471 Real p=0,q=0,r=0,s=0,z=0,t,w,x,y;
00472
00473
00474
00475 Real norm = 0.0;
00476 for (int i = 0; i < nn; i++) {
00477 if ((i < low) || (i > high)) {
00478 d[i] = H[i][i];
00479 e[i] = 0.0;
00480 }
00481 for (int j = max(i-1,0); j < nn; j++) {
00482 norm = norm + abs(H[i][j]);
00483 }
00484 }
00485
00486
00487
00488 int iter = 0;
00489 while (n >= low) {
00490
00491
00492
00493 int l = n;
00494 while (l > low) {
00495 s = abs(H[l-1][l-1]) + abs(H[l][l]);
00496 if (s == 0.0) {
00497 s = norm;
00498 }
00499 if (abs(H[l][l-1]) < eps * s) {
00500 break;
00501 }
00502 l--;
00503 }
00504
00505
00506
00507
00508 if (l == n) {
00509 H[n][n] = H[n][n] + exshift;
00510 d[n] = H[n][n];
00511 e[n] = 0.0;
00512 n--;
00513 iter = 0;
00514
00515
00516
00517 } else if (l == n-1) {
00518 w = H[n][n-1] * H[n-1][n];
00519 p = (H[n-1][n-1] - H[n][n]) / 2.0;
00520 q = p * p + w;
00521 z = sqrt(abs(q));
00522 H[n][n] = H[n][n] + exshift;
00523 H[n-1][n-1] = H[n-1][n-1] + exshift;
00524 x = H[n][n];
00525
00526
00527
00528 if (q >= 0) {
00529 if (p >= 0) {
00530 z = p + z;
00531 } else {
00532 z = p - z;
00533 }
00534 d[n-1] = x + z;
00535 d[n] = d[n-1];
00536 if (z != 0.0) {
00537 d[n] = x - w / z;
00538 }
00539 e[n-1] = 0.0;
00540 e[n] = 0.0;
00541 x = H[n][n-1];
00542 s = abs(x) + abs(z);
00543 p = x / s;
00544 q = z / s;
00545 r = sqrt(p * p+q * q);
00546 p = p / r;
00547 q = q / r;
00548
00549
00550
00551 for (int j = n-1; j < nn; j++) {
00552 z = H[n-1][j];
00553 H[n-1][j] = q * z + p * H[n][j];
00554 H[n][j] = q * H[n][j] - p * z;
00555 }
00556
00557
00558
00559 for (int i = 0; i <= n; i++) {
00560 z = H[i][n-1];
00561 H[i][n-1] = q * z + p * H[i][n];
00562 H[i][n] = q * H[i][n] - p * z;
00563 }
00564
00565
00566
00567 for (int i = low; i <= high; i++) {
00568 z = V[i][n-1];
00569 V[i][n-1] = q * z + p * V[i][n];
00570 V[i][n] = q * V[i][n] - p * z;
00571 }
00572
00573
00574
00575 } else {
00576 d[n-1] = x + p;
00577 d[n] = x + p;
00578 e[n-1] = z;
00579 e[n] = -z;
00580 }
00581 n = n - 2;
00582 iter = 0;
00583
00584
00585
00586 } else {
00587
00588
00589
00590 x = H[n][n];
00591 y = 0.0;
00592 w = 0.0;
00593 if (l < n) {
00594 y = H[n-1][n-1];
00595 w = H[n][n-1] * H[n-1][n];
00596 }
00597
00598
00599
00600 if (iter == 10) {
00601 exshift += x;
00602 for (int i = low; i <= n; i++) {
00603 H[i][i] -= x;
00604 }
00605 s = abs(H[n][n-1]) + abs(H[n-1][n-2]);
00606 x = y = 0.75 * s;
00607 w = -0.4375 * s * s;
00608 }
00609
00610
00611
00612 if (iter == 30) {
00613 s = (y - x) / 2.0;
00614 s = s * s + w;
00615 if (s > 0) {
00616 s = sqrt(s);
00617 if (y < x) {
00618 s = -s;
00619 }
00620 s = x - w / ((y - x) / 2.0 + s);
00621 for (int i = low; i <= n; i++) {
00622 H[i][i] -= s;
00623 }
00624 exshift += s;
00625 x = y = w = 0.964;
00626 }
00627 }
00628
00629 iter = iter + 1;
00630
00631
00632
00633 int m = n-2;
00634 while (m >= l) {
00635 z = H[m][m];
00636 r = x - z;
00637 s = y - z;
00638 p = (r * s - w) / H[m+1][m] + H[m][m+1];
00639 q = H[m+1][m+1] - z - r - s;
00640 r = H[m+2][m+1];
00641 s = abs(p) + abs(q) + abs(r);
00642 p = p / s;
00643 q = q / s;
00644 r = r / s;
00645 if (m == l) {
00646 break;
00647 }
00648 if (abs(H[m][m-1]) * (abs(q) + abs(r)) <
00649 eps * (abs(p) * (abs(H[m-1][m-1]) + abs(z) +
00650 abs(H[m+1][m+1])))) {
00651 break;
00652 }
00653 m--;
00654 }
00655
00656 for (int i = m+2; i <= n; i++) {
00657 H[i][i-2] = 0.0;
00658 if (i > m+2) {
00659 H[i][i-3] = 0.0;
00660 }
00661 }
00662
00663
00664
00665 for (int k = m; k <= n-1; k++) {
00666 int notlast = (k != n-1);
00667 if (k != m) {
00668 p = H[k][k-1];
00669 q = H[k+1][k-1];
00670 r = (notlast ? H[k+2][k-1] : 0.0);
00671 x = abs(p) + abs(q) + abs(r);
00672 if (x != 0.0) {
00673 p = p / x;
00674 q = q / x;
00675 r = r / x;
00676 }
00677 }
00678 if (x == 0.0) {
00679 break;
00680 }
00681 s = sqrt(p * p + q * q + r * r);
00682 if (p < 0) {
00683 s = -s;
00684 }
00685 if (s != 0) {
00686 if (k != m) {
00687 H[k][k-1] = -s * x;
00688 } else if (l != m) {
00689 H[k][k-1] = -H[k][k-1];
00690 }
00691 p = p + s;
00692 x = p / s;
00693 y = q / s;
00694 z = r / s;
00695 q = q / p;
00696 r = r / p;
00697
00698
00699
00700 for (int j = k; j < nn; j++) {
00701 p = H[k][j] + q * H[k+1][j];
00702 if (notlast) {
00703 p = p + r * H[k+2][j];
00704 H[k+2][j] = H[k+2][j] - p * z;
00705 }
00706 H[k][j] = H[k][j] - p * x;
00707 H[k+1][j] = H[k+1][j] - p * y;
00708 }
00709
00710
00711
00712 for (int i = 0; i <= min(n,k+3); i++) {
00713 p = x * H[i][k] + y * H[i][k+1];
00714 if (notlast) {
00715 p = p + z * H[i][k+2];
00716 H[i][k+2] = H[i][k+2] - p * r;
00717 }
00718 H[i][k] = H[i][k] - p;
00719 H[i][k+1] = H[i][k+1] - p * q;
00720 }
00721
00722
00723
00724 for (int i = low; i <= high; i++) {
00725 p = x * V[i][k] + y * V[i][k+1];
00726 if (notlast) {
00727 p = p + z * V[i][k+2];
00728 V[i][k+2] = V[i][k+2] - p * r;
00729 }
00730 V[i][k] = V[i][k] - p;
00731 V[i][k+1] = V[i][k+1] - p * q;
00732 }
00733 }
00734 }
00735 }
00736 }
00737
00738
00739
00740 if (norm == 0.0) {
00741 return;
00742 }
00743
00744 for (n = nn-1; n >= 0; n--) {
00745 p = d[n];
00746 q = e[n];
00747
00748
00749
00750 if (q == 0) {
00751 int l = n;
00752 H[n][n] = 1.0;
00753 for (int i = n-1; i >= 0; i--) {
00754 w = H[i][i] - p;
00755 r = 0.0;
00756 for (int j = l; j <= n; j++) {
00757 r = r + H[i][j] * H[j][n];
00758 }
00759 if (e[i] < 0.0) {
00760 z = w;
00761 s = r;
00762 } else {
00763 l = i;
00764 if (e[i] == 0.0) {
00765 if (w != 0.0) {
00766 H[i][n] = -r / w;
00767 } else {
00768 H[i][n] = -r / (eps * norm);
00769 }
00770
00771
00772
00773 } else {
00774 x = H[i][i+1];
00775 y = H[i+1][i];
00776 q = (d[i] - p) * (d[i] - p) + e[i] * e[i];
00777 t = (x * s - z * r) / q;
00778 H[i][n] = t;
00779 if (abs(x) > abs(z)) {
00780 H[i+1][n] = (-r - w * t) / x;
00781 } else {
00782 H[i+1][n] = (-s - y * t) / z;
00783 }
00784 }
00785
00786
00787
00788 t = abs(H[i][n]);
00789 if ((eps * t) * t > 1) {
00790 for (int j = i; j <= n; j++) {
00791 H[j][n] = H[j][n] / t;
00792 }
00793 }
00794 }
00795 }
00796
00797
00798
00799 } else if (q < 0) {
00800 int l = n-1;
00801
00802
00803
00804 if (abs(H[n][n-1]) > abs(H[n-1][n])) {
00805 H[n-1][n-1] = q / H[n][n-1];
00806 H[n-1][n] = -(H[n][n] - p) / H[n][n-1];
00807 } else {
00808 cdiv(0.0,-H[n-1][n],H[n-1][n-1]-p,q);
00809 H[n-1][n-1] = cdivr;
00810 H[n-1][n] = cdivi;
00811 }
00812 H[n][n-1] = 0.0;
00813 H[n][n] = 1.0;
00814 for (int i = n-2; i >= 0; i--) {
00815 Real ra,sa,vr,vi;
00816 ra = 0.0;
00817 sa = 0.0;
00818 for (int j = l; j <= n; j++) {
00819 ra = ra + H[i][j] * H[j][n-1];
00820 sa = sa + H[i][j] * H[j][n];
00821 }
00822 w = H[i][i] - p;
00823
00824 if (e[i] < 0.0) {
00825 z = w;
00826 r = ra;
00827 s = sa;
00828 } else {
00829 l = i;
00830 if (e[i] == 0) {
00831 cdiv(-ra,-sa,w,q);
00832 H[i][n-1] = cdivr;
00833 H[i][n] = cdivi;
00834 } else {
00835
00836
00837
00838 x = H[i][i+1];
00839 y = H[i+1][i];
00840 vr = (d[i] - p) * (d[i] - p) + e[i] * e[i] - q * q;
00841 vi = (d[i] - p) * 2.0 * q;
00842 if ((vr == 0.0) && (vi == 0.0)) {
00843 vr = eps * norm * (abs(w) + abs(q) +
00844 abs(x) + abs(y) + abs(z));
00845 }
00846 cdiv(x*r-z*ra+q*sa,x*s-z*sa-q*ra,vr,vi);
00847 H[i][n-1] = cdivr;
00848 H[i][n] = cdivi;
00849 if (abs(x) > (abs(z) + abs(q))) {
00850 H[i+1][n-1] = (-ra - w * H[i][n-1] + q * H[i][n]) / x;
00851 H[i+1][n] = (-sa - w * H[i][n] - q * H[i][n-1]) / x;
00852 } else {
00853 cdiv(-r-y*H[i][n-1],-s-y*H[i][n],z,q);
00854 H[i+1][n-1] = cdivr;
00855 H[i+1][n] = cdivi;
00856 }
00857 }
00858
00859
00860
00861 t = max(abs(H[i][n-1]),abs(H[i][n]));
00862 if ((eps * t) * t > 1) {
00863 for (int j = i; j <= n; j++) {
00864 H[j][n-1] = H[j][n-1] / t;
00865 H[j][n] = H[j][n] / t;
00866 }
00867 }
00868 }
00869 }
00870 }
00871 }
00872
00873
00874
00875 for (int i = 0; i < nn; i++) {
00876 if (i < low || i > high) {
00877 for (int j = i; j < nn; j++) {
00878 V[i][j] = H[i][j];
00879 }
00880 }
00881 }
00882
00883
00884
00885 for (int j = nn-1; j >= low; j--) {
00886 for (int i = low; i <= high; i++) {
00887 z = 0.0;
00888 for (int k = low; k <= min(j,high); k++) {
00889 z = z + V[i][k] * H[k][j];
00890 }
00891 V[i][j] = z;
00892 }
00893 }
00894 }
00895
00896 public:
00897
00898
00903 Eigenvalue(const TNT::Array2D<Real> &A) {
00904 n = A.dim2();
00905 V = Array2D<Real>(n,n);
00906 d = Array1D<Real>(n);
00907 e = Array1D<Real>(n);
00908
00909 issymmetric = 1;
00910 for (int j = 0; (j < n) && issymmetric; j++) {
00911 for (int i = 0; (i < n) && issymmetric; i++) {
00912 issymmetric = (A[i][j] == A[j][i]);
00913 }
00914 }
00915
00916 if (issymmetric) {
00917 for (int i = 0; i < n; i++) {
00918 for (int j = 0; j < n; j++) {
00919 V[i][j] = A[i][j];
00920 }
00921 }
00922
00923
00924 tred2();
00925
00926
00927 tql2();
00928
00929 } else {
00930 H = TNT::Array2D<Real>(n,n);
00931 ort = TNT::Array1D<Real>(n);
00932
00933 for (int j = 0; j < n; j++) {
00934 for (int i = 0; i < n; i++) {
00935 H[i][j] = A[i][j];
00936 }
00937 }
00938
00939
00940 orthes();
00941
00942
00943 hqr2();
00944 }
00945 }
00946
00947
00952 void getV (TNT::Array2D<Real> &V_) {
00953 V_ = V;
00954 return;
00955 }
00956
00961 void getRealEigenvalues (TNT::Array1D<Real> &d_) {
00962 d_ = d;
00963 return ;
00964 }
00965
00971 void getImagEigenvalues (TNT::Array1D<Real> &e_) {
00972 e_ = e;
00973 return;
00974 }
00975
00976
01010 void getD (TNT::Array2D<Real> &D) {
01011 D = Array2D<Real>(n,n);
01012 for (int i = 0; i < n; i++) {
01013 for (int j = 0; j < n; j++) {
01014 D[i][j] = 0.0;
01015 }
01016 D[i][i] = d[i];
01017 if (e[i] > 0) {
01018 D[i][i+1] = e[i];
01019 } else if (e[i] < 0) {
01020 D[i][i-1] = e[i];
01021 }
01022 }
01023 }
01024 };
01025
01026 }
01027
01028
01029 #endif
01030